Geant4 10.7.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 void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *)
 
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 void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
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 GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, 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 *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForLossfParticleChange
 
const G4ParticleDefinitionfParticle
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

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 62 of file G4PenelopeBremsstrahlungModel.hh.

Constructor & Destructor Documentation

◆ G4PenelopeBremsstrahlungModel()

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

Definition at line 63 of file G4PenelopeBremsstrahlungModel.cc.

66 isInitialised(false),energyGrid(0),
67 XSTableElectron(0),XSTablePositron(0),fPenelopeFSHelper(0),
68 fPenelopeAngular(0),fLocalTable(false)
69
70{
71 fIntrinsicLowEnergyLimit = 100.0*eV;
72 fIntrinsicHighEnergyLimit = 100.0*GeV;
73 nBins = 200;
74
75 if (part)
76 SetParticle(part);
77
78 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
79 //
81 //
82 verboseLevel= 0;
83
84 // Verbosity scale:
85 // 0 = nothing
86 // 1 = warning for energy non-conservation
87 // 2 = details of energy budget
88 // 3 = calculation of cross sections, file openings, sampling of atoms
89 // 4 = entering in methods
90
91 // Atomic deexcitation model activated by default
93
94}
static G4PenelopeOscillatorManager * GetOscillatorManager()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813

◆ ~G4PenelopeBremsstrahlungModel()

G4PenelopeBremsstrahlungModel::~G4PenelopeBremsstrahlungModel ( )
virtual

Definition at line 98 of file G4PenelopeBremsstrahlungModel.cc.

99{
100 if (IsMaster() || fLocalTable)
101 {
102 ClearTables();
103 if (fPenelopeFSHelper)
104 delete fPenelopeFSHelper;
105 }
106 // This is thread-local at the moment
107 if (fPenelopeAngular)
108 delete fPenelopeAngular;
109}
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

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 281 of file G4PenelopeBremsstrahlungModel.cc.

287{
288 G4cout << "*** G4PenelopeBremsstrahlungModel -- WARNING ***" << G4endl;
289 G4cout << "Penelope Bremsstrahlung model v2008 does not calculate cross section _per atom_ " << G4endl;
290 G4cout << "so the result is always zero. For physics values, please invoke " << G4endl;
291 G4cout << "GetCrossSectionPerVolume() or GetMeanFreePath() via the G4EmCalculator" << G4endl;
292 return 0;
293}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 297 of file G4PenelopeBremsstrahlungModel.cc.

301{
302 if (verboseLevel > 3)
303 G4cout << "Calling ComputeDEDX() of G4PenelopeBremsstrahlungModel" << G4endl;
304
305 G4PenelopeCrossSection* theXS = GetCrossSectionTableForCouple(theParticle,material,
306 cutEnergy);
307 G4double sPowerPerMolecule = 0.0;
308 if (theXS)
309 sPowerPerMolecule = theXS->GetSoftStoppingPower(kineticEnergy);
310
311
312 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
313 G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
314
315 G4double moleculeDensity = 0.;
316 if (atPerMol)
317 moleculeDensity = atomDensity/atPerMol;
318
319 G4double sPowerPerVolume = sPowerPerMolecule*moleculeDensity;
320
321 if (verboseLevel > 2)
322 {
323 G4cout << "G4PenelopeBremsstrahlungModel " << G4endl;
324 G4cout << "Stopping power < " << cutEnergy/keV << " keV at " <<
325 kineticEnergy/keV << " keV = " <<
326 sPowerPerVolume/(keV/mm) << " keV/mm" << G4endl;
327 }
328 return sPowerPerVolume;
329}
double G4double
Definition: G4Types.hh:83
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:207
G4double GetSoftStoppingPower(G4double energy) const
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 233 of file G4PenelopeBremsstrahlungModel.cc.

238{
239 //
240 if (verboseLevel > 3)
241 G4cout << "Calling CrossSectionPerVolume() of G4PenelopeBremsstrahlungModel" << G4endl;
242
243 SetupForMaterial(theParticle, material, energy);
244
245 G4double crossPerMolecule = 0.;
246
247 G4PenelopeCrossSection* theXS = GetCrossSectionTableForCouple(theParticle,material,
248 cutEnergy);
249
250 if (theXS)
251 crossPerMolecule = theXS->GetHardCrossSection(energy);
252
253 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
254 G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
255
256 if (verboseLevel > 3)
257 G4cout << "Material " << material->GetName() << " has " << atPerMol <<
258 "atoms per molecule" << G4endl;
259
260 G4double moleculeDensity = 0.;
261 if (atPerMol)
262 moleculeDensity = atomDensity/atPerMol;
263
264 G4double crossPerVolume = crossPerMolecule*moleculeDensity;
265
266 if (verboseLevel > 2)
267 {
268 G4cout << "G4PenelopeBremsstrahlungModel " << G4endl;
269 G4cout << "Mean free path for gamma emission > " << cutEnergy/keV << " keV at " <<
270 energy/keV << " keV = " << (1./crossPerVolume)/mm << " mm" << G4endl;
271 }
272
273 return crossPerVolume;
274}
const G4String & GetName() const
Definition: G4Material.hh:175
G4double GetHardCrossSection(G4double energy) const
Returns hard cross section at the given energy.
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:449
G4double energy(const ThreeVector &p, const G4double m)

◆ GetVerbosityLevel()

G4int G4PenelopeBremsstrahlungModel::GetVerbosityLevel ( )
inline

Definition at line 108 of file G4PenelopeBremsstrahlungModel.hh.

108{return verboseLevel;};

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 113 of file G4PenelopeBremsstrahlungModel.cc.

115{
116 if (verboseLevel > 3)
117 G4cout << "Calling G4PenelopeBremsstrahlungModel::Initialise()" << G4endl;
118
119 SetParticle(part);
120
121 if (IsMaster() && part == fParticle)
122 {
123
124 if (!fPenelopeFSHelper)
125 fPenelopeFSHelper = new G4PenelopeBremsstrahlungFS(verboseLevel);
126 if (!fPenelopeAngular)
127 fPenelopeAngular = new G4PenelopeBremsstrahlungAngular();
128 //Clear and re-build the tables
129 ClearTables();
130
131 //forces the cleaning of tables, in this specific case
132 if (fPenelopeAngular)
133 fPenelopeAngular->Initialize();
134
135 //Set the number of bins for the tables. 20 points per decade
136 nBins = (size_t) (20*std::log10(HighEnergyLimit()/LowEnergyLimit()));
137 nBins = std::max(nBins,(size_t)100);
138 energyGrid = new G4PhysicsLogVector(LowEnergyLimit(),
140 nBins-1); //one hidden bin is added
141
142
143 XSTableElectron = new
144 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*>;
145 XSTablePositron = new
146 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*>;
147
148 G4ProductionCutsTable* theCoupleTable =
150
151 //Build tables for all materials
152 for (size_t i=0;i<theCoupleTable->GetTableSize();i++)
153 {
154 const G4Material* theMat =
155 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
156 //Forces the building of the helper tables
157 fPenelopeFSHelper->BuildScaledXSTable(theMat,theCuts.at(i),IsMaster());
158 fPenelopeAngular->PrepareTables(theMat,IsMaster());
159 BuildXSTable(theMat,theCuts.at(i));
160
161 }
162
163
164 if (verboseLevel > 2) {
165 G4cout << "Penelope Bremsstrahlung model v2008 is initialized " << G4endl
166 << "Energy range: "
167 << LowEnergyLimit() / keV << " keV - "
168 << HighEnergyLimit() / GeV << " GeV."
169 << G4endl;
170 }
171 }
172
173 if(isInitialised) return;
175 isInitialised = true;
176}
const G4Material * GetMaterial() const
void PrepareTables(const G4Material *material, G4bool isMaster)
Reserved for Master Model.
void BuildScaledXSTable(const G4Material *material, G4double cut, G4bool isMaster)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118

◆ InitialiseLocal()

void G4PenelopeBremsstrahlungModel::InitialiseLocal ( const G4ParticleDefinition part,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 179 of file G4PenelopeBremsstrahlungModel.cc.

181{
182 if (verboseLevel > 3)
183 G4cout << "Calling G4PenelopeBremsstrahlungModel::InitialiseLocal()" << G4endl;
184
185 //
186 //Check that particle matches: one might have multiple master models (e.g.
187 //for e+ and e-).
188 //
189 if (part == fParticle)
190 {
191 //Get the const table pointers from the master to the workers
192 const G4PenelopeBremsstrahlungModel* theModel =
193 static_cast<G4PenelopeBremsstrahlungModel*> (masterModel);
194
195 //Copy pointers to the data tables
196 energyGrid = theModel->energyGrid;
197 XSTableElectron = theModel->XSTableElectron;
198 XSTablePositron = theModel->XSTablePositron;
199 fPenelopeFSHelper = theModel->fPenelopeFSHelper;
200
201 //fPenelopeAngular = theModel->fPenelopeAngular;
202
203 //created in each thread and initialized.
204 if (!fPenelopeAngular)
205 fPenelopeAngular = new G4PenelopeBremsstrahlungAngular();
206 //forces the cleaning of tables, in this specific case
207 if (fPenelopeAngular)
208 fPenelopeAngular->Initialize();
209
210 G4ProductionCutsTable* theCoupleTable =
212 //Build tables for all materials
213 for (size_t i=0;i<theCoupleTable->GetTableSize();i++)
214 {
215 const G4Material* theMat =
216 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
217 fPenelopeAngular->PrepareTables(theMat,IsMaster());
218 }
219
220
221 //copy the data
222 nBins = theModel->nBins;
223
224 //Same verbosity for all workers, as the master
225 verboseLevel = theModel->verboseLevel;
226 }
227
228 return;
229}

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 472 of file G4PenelopeBremsstrahlungModel.cc.

474{
475 return fIntrinsicLowEnergyLimit;
476}

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 333 of file G4PenelopeBremsstrahlungModel.cc.

338{
339 if (verboseLevel > 3)
340 G4cout << "Calling SampleSecondaries() of G4PenelopeBremsstrahlungModel" << G4endl;
341
342 G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
343 const G4Material* material = couple->GetMaterial();
344
345 if (kineticEnergy <= fIntrinsicLowEnergyLimit)
346 {
349 return ;
350 }
351
352 G4ParticleMomentum particleDirection0 = aDynamicParticle->GetMomentumDirection();
353 //This is the momentum
354 G4ThreeVector initialMomentum = aDynamicParticle->GetMomentum();
355
356 //Not enough energy to produce a secondary! Return with nothing happened
357 if (kineticEnergy < cutG)
358 return;
359
360 if (verboseLevel > 3)
361 G4cout << "Going to sample gamma energy for: " <<material->GetName() << " " <<
362 "energy = " << kineticEnergy/keV << ", cut = " << cutG/keV << G4endl;
363
364 //Sample gamma's energy according to the spectrum
365 G4double gammaEnergy =
366 fPenelopeFSHelper->SampleGammaEnergy(kineticEnergy,material,cutG);
367
368 if (verboseLevel > 3)
369 G4cout << "Sampled gamma energy: " << gammaEnergy/keV << " keV" << G4endl;
370
371 //Now sample the direction for the Gamma. Notice that the rotation is already done
372 //Z is unused here, I plug 0. The information is in the material pointer
373 G4ThreeVector gammaDirection1 =
374 fPenelopeAngular->SampleDirection(aDynamicParticle,gammaEnergy,0,material);
375
376 if (verboseLevel > 3)
377 G4cout << "Sampled cosTheta for e-: " << gammaDirection1.cosTheta() << G4endl;
378
379 G4double residualPrimaryEnergy = kineticEnergy-gammaEnergy;
380 if (residualPrimaryEnergy < 0)
381 {
382 //Ok we have a problem, all energy goes with the gamma
383 gammaEnergy += residualPrimaryEnergy;
384 residualPrimaryEnergy = 0.0;
385 }
386
387 //Produce final state according to momentum conservation
388 G4ThreeVector particleDirection1 = initialMomentum - gammaEnergy*gammaDirection1;
389 particleDirection1 = particleDirection1.unit(); //normalize
390
391 //Update the primary particle
392 if (residualPrimaryEnergy > 0.)
393 {
394 fParticleChange->ProposeMomentumDirection(particleDirection1);
395 fParticleChange->SetProposedKineticEnergy(residualPrimaryEnergy);
396 }
397 else
398 {
400 }
401
402 //Now produce the photon
404 gammaDirection1,
405 gammaEnergy);
406 fvect->push_back(theGamma);
407
408 if (verboseLevel > 1)
409 {
410 G4cout << "-----------------------------------------------------------" << G4endl;
411 G4cout << "Energy balance from G4PenelopeBremsstrahlung" << G4endl;
412 G4cout << "Incoming primary energy: " << kineticEnergy/keV << " keV" << G4endl;
413 G4cout << "-----------------------------------------------------------" << G4endl;
414 G4cout << "Outgoing primary energy: " << residualPrimaryEnergy/keV << " keV" << G4endl;
415 G4cout << "Bremsstrahlung photon " << gammaEnergy/keV << " keV" << G4endl;
416 G4cout << "Total final state: " << (residualPrimaryEnergy+gammaEnergy)/keV
417 << " keV" << G4endl;
418 G4cout << "-----------------------------------------------------------" << G4endl;
419 }
420
421 if (verboseLevel > 0)
422 {
423 G4double energyDiff = std::fabs(residualPrimaryEnergy+gammaEnergy-kineticEnergy);
424 if (energyDiff > 0.05*keV)
425 G4cout << "Warning from G4PenelopeBremsstrahlung: problem with energy conservation: "
426 <<
427 (residualPrimaryEnergy+gammaEnergy)/keV <<
428 " keV (final) vs. " <<
429 kineticEnergy/keV << " keV (initial)" << G4endl;
430 }
431 return;
432}
Hep3Vector unit() const
double cosTheta() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
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)
Samples the direction of the outgoing photon (in global coordinates).
G4double SampleGammaEnergy(G4double energy, const G4Material *, const G4double cut) const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetVerbosityLevel()

void G4PenelopeBremsstrahlungModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 107 of file G4PenelopeBremsstrahlungModel.hh.

107{verboseLevel = lev;};

Member Data Documentation

◆ fParticle

const G4ParticleDefinition* G4PenelopeBremsstrahlungModel::fParticle
protected

Definition at line 112 of file G4PenelopeBremsstrahlungModel.hh.

Referenced by Initialise(), and InitialiseLocal().

◆ fParticleChange

G4ParticleChangeForLoss* G4PenelopeBremsstrahlungModel::fParticleChange
protected

Definition at line 111 of file G4PenelopeBremsstrahlungModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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