Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeBremsstrahlungModel Class Reference

#include <G4PenelopeBremsstrahlungModel.hh>

+ Inheritance diagram for G4PenelopeBremsstrahlungModel:

Public Member Functions

 G4PenelopeBremsstrahlungModel (const G4ParticleDefinition *p=0, const G4String &processName="PenBrem")
 
virtual ~G4PenelopeBremsstrahlungModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *theParticle, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *theParticle, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
void SetVerbosityLevel (G4int lev)
 
G4int GetVerbosityLevel ()
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForLossfParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 63 of file G4PenelopeBremsstrahlungModel.hh.

Constructor & Destructor Documentation

◆ G4PenelopeBremsstrahlungModel()

G4PenelopeBremsstrahlungModel::G4PenelopeBremsstrahlungModel ( const G4ParticleDefinition p = 0,
const G4String processName = "PenBrem" 
)

Definition at line 59 of file G4PenelopeBremsstrahlungModel.cc.

61 :G4VEmModel(nam),fParticleChange(0),isInitialised(false),energyGrid(0),
62 XSTableElectron(0),XSTablePositron(0)
63
64{
65 fIntrinsicLowEnergyLimit = 100.0*eV;
66 fIntrinsicHighEnergyLimit = 100.0*GeV;
67 nBins = 200;
68
69 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
70 //
72 //
73 verboseLevel= 0;
74
75 // Verbosity scale:
76 // 0 = nothing
77 // 1 = warning for energy non-conservation
78 // 2 = details of energy budget
79 // 3 = calculation of cross sections, file openings, sampling of atoms
80 // 4 = entering in methods
81
82 fPenelopeFSHelper = new G4PenelopeBremsstrahlungFS();
83 fPenelopeAngular = new G4PenelopeBremsstrahlungAngular();
84
85 // Atomic deexcitation model activated by default
87
88}
static G4PenelopeOscillatorManager * GetOscillatorManager()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641

◆ ~G4PenelopeBremsstrahlungModel()

G4PenelopeBremsstrahlungModel::~G4PenelopeBremsstrahlungModel ( )
virtual

Definition at line 92 of file G4PenelopeBremsstrahlungModel.cc.

93{
94 ClearTables();
95 if (fPenelopeFSHelper)
96 delete fPenelopeFSHelper;
97 if (fPenelopeAngular)
98 delete fPenelopeAngular;
99}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4PenelopeBremsstrahlungModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition theParticle,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 193 of file G4PenelopeBremsstrahlungModel.cc.

199{
200 G4cout << "*** G4PenelopeBremsstrahlungModel -- WARNING ***" << G4endl;
201 G4cout << "Penelope Bremsstrahlung model v2008 does not calculate cross section _per atom_ " << G4endl;
202 G4cout << "so the result is always zero. For physics values, please invoke " << G4endl;
203 G4cout << "GetCrossSectionPerVolume() or GetMeanFreePath() via the G4EmCalculator" << G4endl;
204 return 0;
205}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ ComputeDEDXPerVolume()

G4double G4PenelopeBremsstrahlungModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition theParticle,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 209 of file G4PenelopeBremsstrahlungModel.cc.

213{
214 if (verboseLevel > 3)
215 G4cout << "Calling ComputeDEDX() of G4PenelopeBremsstrahlungModel" << G4endl;
216
217 G4PenelopeCrossSection* theXS = GetCrossSectionTableForCouple(theParticle,material,
218 cutEnergy);
219 G4double sPowerPerMolecule = 0.0;
220 if (theXS)
221 sPowerPerMolecule = theXS->GetSoftStoppingPower(kineticEnergy);
222
223
224 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
225 G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
226
227 G4double moleculeDensity = 0.;
228 if (atPerMol)
229 moleculeDensity = atomDensity/atPerMol;
230
231 G4double sPowerPerVolume = sPowerPerMolecule*moleculeDensity;
232
233 if (verboseLevel > 2)
234 {
235 G4cout << "G4PenelopeBremsstrahlungModel " << G4endl;
236 G4cout << "Stopping power < " << cutEnergy/keV << " keV at " <<
237 kineticEnergy/keV << " keV = " <<
238 sPowerPerVolume/(keV/mm) << " keV/mm" << G4endl;
239 }
240 return sPowerPerVolume;
241}
double G4double
Definition: G4Types.hh:64
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:208
G4double GetSoftStoppingPower(G4double energy)
Returns the total stopping power due to soft collisions.
G4double GetAtomsPerMolecule(const G4Material *)
Returns the total number of atoms per molecule.

◆ CrossSectionPerVolume()

G4double G4PenelopeBremsstrahlungModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition theParticle,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 145 of file G4PenelopeBremsstrahlungModel.cc.

150{
151 //
152 if (verboseLevel > 3)
153 G4cout << "Calling CrossSectionPerVolume() of G4PenelopeBremsstrahlungModel" << G4endl;
154
155 SetupForMaterial(theParticle, material, energy);
156
157 G4double crossPerMolecule = 0.;
158
159 G4PenelopeCrossSection* theXS = GetCrossSectionTableForCouple(theParticle,material,
160 cutEnergy);
161
162 if (theXS)
163 crossPerMolecule = theXS->GetHardCrossSection(energy);
164
165 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
166 G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
167
168 if (verboseLevel > 3)
169 G4cout << "Material " << material->GetName() << " has " << atPerMol <<
170 "atoms per molecule" << G4endl;
171
172 G4double moleculeDensity = 0.;
173 if (atPerMol)
174 moleculeDensity = atomDensity/atPerMol;
175
176 G4double crossPerVolume = crossPerMolecule*moleculeDensity;
177
178 if (verboseLevel > 2)
179 {
180 G4cout << "G4PenelopeBremsstrahlungModel " << G4endl;
181 G4cout << "Mean free path for gamma emission > " << cutEnergy/keV << " keV at " <<
182 energy/keV << " keV = " << (1./crossPerVolume)/mm << " mm" << G4endl;
183 }
184
185 return crossPerVolume;
186}
const G4String & GetName() const
Definition: G4Material.hh:177
G4double GetHardCrossSection(G4double energy)
Returns hard cross section at the given energy.
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:311

◆ GetVerbosityLevel()

G4int G4PenelopeBremsstrahlungModel::GetVerbosityLevel ( )
inline

Definition at line 107 of file G4PenelopeBremsstrahlungModel.hh.

107{return verboseLevel;};

◆ Initialise()

void G4PenelopeBremsstrahlungModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 103 of file G4PenelopeBremsstrahlungModel.cc.

105{
106 if (verboseLevel > 3)
107 G4cout << "Calling G4PenelopeBremsstrahlungModel::Initialise()" << G4endl;
108
109 //Clear and re-build the tables
110 ClearTables();
111
112 //forces the cleaning of tables, in this specific case
113 if (fPenelopeAngular)
114 fPenelopeAngular->Initialize();
115
116
117 //Set the number of bins for the tables. 20 points per decade
118 nBins = (size_t) (20*std::log10(HighEnergyLimit()/LowEnergyLimit()));
119 nBins = std::max(nBins,(size_t)100);
120 energyGrid = new G4PhysicsLogVector(LowEnergyLimit(),
122 nBins-1); //one hidden bin is added
123
124
125 XSTableElectron = new
126 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*>;
127 XSTablePositron = new
128 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*>;
129
130 if (verboseLevel > 2) {
131 G4cout << "Penelope Bremsstrahlung model v2008 is initialized " << G4endl
132 << "Energy range: "
133 << LowEnergyLimit() / keV << " keV - "
134 << HighEnergyLimit() / GeV << " GeV."
135 << G4endl;
136 }
137
138 if(isInitialised) return;
140 isInitialised = true;
141}
void Initialize()
The Initialize() method forces the cleaning of tables.
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95

◆ MinEnergyCut()

G4double G4PenelopeBremsstrahlungModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple  
)
virtual

Definition at line 386 of file G4PenelopeBremsstrahlungModel.cc.

388{
389 return fIntrinsicLowEnergyLimit;
390}

◆ SampleSecondaries()

void G4PenelopeBremsstrahlungModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicParticle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 245 of file G4PenelopeBremsstrahlungModel.cc.

250{
251 if (verboseLevel > 3)
252 G4cout << "Calling SampleSecondaries() of G4PenelopeBremsstrahlungModel" << G4endl;
253
254 G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
255 const G4Material* material = couple->GetMaterial();
256
257 if (kineticEnergy <= fIntrinsicLowEnergyLimit)
258 {
261 return ;
262 }
263
264 G4ParticleMomentum particleDirection0 = aDynamicParticle->GetMomentumDirection();
265 //This is the momentum
266 G4ThreeVector initialMomentum = aDynamicParticle->GetMomentum();
267
268 //Not enough energy to produce a secondary! Return with nothing happened
269 if (kineticEnergy < cutG)
270 return;
271
272 if (verboseLevel > 3)
273 G4cout << "Going to sample gamma energy for: " <<material->GetName() << " " <<
274 "energy = " << kineticEnergy/keV << ", cut = " << cutG/keV << G4endl;
275
276 //Sample gamma's energy according to the spectrum
277 G4double gammaEnergy =
278 fPenelopeFSHelper->SampleGammaEnergy(kineticEnergy,material,cutG);
279
280 if (verboseLevel > 3)
281 G4cout << "Sampled gamma energy: " << gammaEnergy/keV << " keV" << G4endl;
282
283 //Now sample the direction for the Gamma. Notice that the rotation is already done
284 //Z is unused here, I plug 0. The information is in the material pointer
285 G4ThreeVector gammaDirection1 =
286 fPenelopeAngular->SampleDirection(aDynamicParticle,gammaEnergy,0,material);
287
288 if (verboseLevel > 3)
289 G4cout << "Sampled cosTheta for e-: " << gammaDirection1.cosTheta() << G4endl;
290
291 G4double residualPrimaryEnergy = kineticEnergy-gammaEnergy;
292 if (residualPrimaryEnergy < 0)
293 {
294 //Ok we have a problem, all energy goes with the gamma
295 gammaEnergy += residualPrimaryEnergy;
296 residualPrimaryEnergy = 0.0;
297 }
298
299 //Produce final state according to momentum conservation
300 G4ThreeVector particleDirection1 = initialMomentum - gammaEnergy*gammaDirection1;
301 particleDirection1 = particleDirection1.unit(); //normalize
302
303 //Update the primary particle
304 if (residualPrimaryEnergy > 0.)
305 {
306 fParticleChange->ProposeMomentumDirection(particleDirection1);
307 fParticleChange->SetProposedKineticEnergy(residualPrimaryEnergy);
308 }
309 else
310 {
312 }
313
314 //Now produce the photon
316 gammaDirection1,
317 gammaEnergy);
318 fvect->push_back(theGamma);
319
320 if (verboseLevel > 1)
321 {
322 G4cout << "-----------------------------------------------------------" << G4endl;
323 G4cout << "Energy balance from G4PenelopeBremsstrahlung" << G4endl;
324 G4cout << "Incoming primary energy: " << kineticEnergy/keV << " keV" << G4endl;
325 G4cout << "-----------------------------------------------------------" << G4endl;
326 G4cout << "Outgoing primary energy: " << residualPrimaryEnergy/keV << " keV" << G4endl;
327 G4cout << "Bremsstrahlung photon " << gammaEnergy/keV << " keV" << G4endl;
328 G4cout << "Total final state: " << (residualPrimaryEnergy+gammaEnergy)/keV
329 << " keV" << G4endl;
330 G4cout << "-----------------------------------------------------------" << G4endl;
331 }
332
333 if (verboseLevel > 0)
334 {
335 G4double energyDiff = std::fabs(residualPrimaryEnergy+gammaEnergy-kineticEnergy);
336 if (energyDiff > 0.05*keV)
337 G4cout << "Warning from G4PenelopeBremsstrahlung: problem with energy conservation: "
338 <<
339 (residualPrimaryEnergy+gammaEnergy)/keV <<
340 " keV (final) vs. " <<
341 kineticEnergy/keV << " keV (initial)" << G4endl;
342 }
343 return;
344}
Hep3Vector unit() const
double cosTheta() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=0)
G4double SampleGammaEnergy(G4double energy, const G4Material *, G4double cut)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetVerbosityLevel()

void G4PenelopeBremsstrahlungModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 106 of file G4PenelopeBremsstrahlungModel.hh.

106{verboseLevel = lev;};

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForLoss* G4PenelopeBremsstrahlungModel::fParticleChange
protected

Definition at line 110 of file G4PenelopeBremsstrahlungModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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