Geant4 11.1.1
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=nullptr, const G4String &processName="PenBrem")
 
virtual ~G4PenelopeBremsstrahlungModel ()
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *theParticle, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *theParticle, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
 
void SetVerbosityLevel (G4int lev)
 
G4int GetVerbosityLevel ()
 
G4PenelopeBremsstrahlungModeloperator= (const G4PenelopeBremsstrahlungModel &right)=delete
 
 G4PenelopeBremsstrahlungModel (const G4PenelopeBremsstrahlungModel &)=delete
 
- 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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
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)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
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 *)
 
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 = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

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() [1/2]

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

Definition at line 63 of file G4PenelopeBremsstrahlungModel.cc.

65 :G4VEmModel(nam),fParticleChange(nullptr),fParticle(nullptr),
66 fPenelopeFSHelper(nullptr),fPenelopeAngular(nullptr),fEnergyGrid(nullptr),
67 fXSTableElectron(nullptr),fXSTablePositron(nullptr),
68 fIsInitialised(false),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 fVerboseLevel= 0;
83 // Verbosity scale:
84 // 0 = nothing
85 // 1 = warning for energy non-conservation
86 // 2 = details of energy budget
87 // 3 = calculation of cross sections, file openings, sampling of atoms
88 // 4 = entering in methods
89
90 // Atomic deexcitation model activated by default
92
93}
static G4PenelopeOscillatorManager * GetOscillatorManager()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802

◆ ~G4PenelopeBremsstrahlungModel()

G4PenelopeBremsstrahlungModel::~G4PenelopeBremsstrahlungModel ( )
virtual

Definition at line 97 of file G4PenelopeBremsstrahlungModel.cc.

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

◆ G4PenelopeBremsstrahlungModel() [2/2]

G4PenelopeBremsstrahlungModel::G4PenelopeBremsstrahlungModel ( const G4PenelopeBremsstrahlungModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 271 of file G4PenelopeBremsstrahlungModel.cc.

277{
278 G4cout << "*** G4PenelopeBremsstrahlungModel -- WARNING ***" << G4endl;
279 G4cout << "Penelope Bremsstrahlung model v2008 does not calculate cross section _per atom_ " << G4endl;
280 G4cout << "so the result is always zero. For physics values, please invoke " << G4endl;
281 G4cout << "GetCrossSectionPerVolume() or GetMeanFreePath() via the G4EmCalculator" << G4endl;
282 return 0;
283}
#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 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 287 of file G4PenelopeBremsstrahlungModel.cc.

291{
292 if (fVerboseLevel > 3)
293 G4cout << "Calling ComputeDEDX() of G4PenelopeBremsstrahlungModel" << G4endl;
294
295 G4PenelopeCrossSection* theXS = GetCrossSectionTableForCouple(theParticle,material,
296 cutEnergy);
297 G4double sPowerPerMolecule = 0.0;
298 if (theXS)
299 sPowerPerMolecule = theXS->GetSoftStoppingPower(kineticEnergy);
300
301 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
302 G4double atPerMol = fOscManager->GetAtomsPerMolecule(material);
303
304 G4double moleculeDensity = 0.;
305 if (atPerMol)
306 moleculeDensity = atomDensity/atPerMol;
307
308 G4double sPowerPerVolume = sPowerPerMolecule*moleculeDensity;
309
310 if (fVerboseLevel > 2)
311 {
312 G4cout << "G4PenelopeBremsstrahlungModel " << G4endl;
313 G4cout << "Stopping power < " << cutEnergy/keV << " keV at " <<
314 kineticEnergy/keV << " keV = " <<
315 sPowerPerVolume/(keV/mm) << " keV/mm" << G4endl;
316 }
317 return sPowerPerVolume;
318}
double G4double
Definition: G4Types.hh:83
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
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 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 225 of file G4PenelopeBremsstrahlungModel.cc.

230{
231 //
232 if (fVerboseLevel > 3)
233 G4cout << "Calling CrossSectionPerVolume() of G4PenelopeBremsstrahlungModel" << G4endl;
234
235 SetupForMaterial(theParticle, material, energy);
236 G4double crossPerMolecule = 0.;
237
238 G4PenelopeCrossSection* theXS = GetCrossSectionTableForCouple(theParticle,material,
239 cutEnergy);
240 if (theXS)
241 crossPerMolecule = theXS->GetHardCrossSection(energy);
242
243 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
244 G4double atPerMol = fOscManager->GetAtomsPerMolecule(material);
245
246 if (fVerboseLevel > 3)
247 G4cout << "Material " << material->GetName() << " has " << atPerMol <<
248 "atoms per molecule" << G4endl;
249
250 G4double moleculeDensity = 0.;
251 if (atPerMol)
252 moleculeDensity = atomDensity/atPerMol;
253
254 G4double crossPerVolume = crossPerMolecule*moleculeDensity;
255
256 if (fVerboseLevel > 2)
257 {
258 G4cout << "G4PenelopeBremsstrahlungModel " << G4endl;
259 G4cout << "Mean free path for gamma emission > " << cutEnergy/keV << " keV at " <<
260 energy/keV << " keV = " << (1./crossPerVolume)/mm << " mm" << G4endl;
261 }
262
263 return crossPerVolume;
264}
const G4String & GetName() const
Definition: G4Material.hh:172
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:383
G4double energy(const ThreeVector &p, const G4double m)

◆ GetVerbosityLevel()

G4int G4PenelopeBremsstrahlungModel::GetVerbosityLevel ( )
inline

Definition at line 104 of file G4PenelopeBremsstrahlungModel.hh.

104{return fVerboseLevel;};

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 112 of file G4PenelopeBremsstrahlungModel.cc.

114{
115 if (fVerboseLevel > 3)
116 G4cout << "Calling G4PenelopeBremsstrahlungModel::Initialise()" << G4endl;
117
118 SetParticle(part);
119
120 if (IsMaster() && part == fParticle)
121 {
122 if (!fPenelopeFSHelper)
123 fPenelopeFSHelper = new G4PenelopeBremsstrahlungFS(fVerboseLevel);
124 if (!fPenelopeAngular)
125 fPenelopeAngular = new G4PenelopeBremsstrahlungAngular();
126 //Clear and re-build the tables
127 ClearTables();
128
129 //forces the cleaning of tables, in this specific case
130 if (fPenelopeAngular)
131 fPenelopeAngular->Initialize();
132
133 //Set the number of bins for the tables. 20 points per decade
134 nBins = (std::size_t) (20*std::log10(HighEnergyLimit()/LowEnergyLimit()));
135 nBins = std::max(nBins,(std::size_t)100);
136 fEnergyGrid = new G4PhysicsLogVector(LowEnergyLimit(),
138 nBins-1); //one hidden bin is added
139
140 fXSTableElectron = new
141 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*>;
142 fXSTablePositron = new
143 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*>;
144
145 G4ProductionCutsTable* theCoupleTable =
147
148 //Build tables for all materials
149 for (G4int i=0;i<(G4int)theCoupleTable->GetTableSize();++i)
150 {
151 const G4Material* theMat =
152 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
153 //Forces the building of the helper tables
154 fPenelopeFSHelper->BuildScaledXSTable(theMat,theCuts.at(i),IsMaster());
155 fPenelopeAngular->PrepareTables(theMat,IsMaster());
156 BuildXSTable(theMat,theCuts.at(i));
157
158 }
159
160 if (fVerboseLevel > 2) {
161 G4cout << "Penelope Bremsstrahlung model v2008 is initialized " << G4endl
162 << "Energy range: "
163 << LowEnergyLimit() / keV << " keV - "
164 << HighEnergyLimit() / GeV << " GeV."
165 << G4endl;
166 }
167 }
168
169 if(fIsInitialised) return;
171 fIsInitialised = true;
172}
int G4int
Definition: G4Types.hh:85
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:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 176 of file G4PenelopeBremsstrahlungModel.cc.

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

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 461 of file G4PenelopeBremsstrahlungModel.cc.

463{
464 return fIntrinsicLowEnergyLimit;
465}

◆ operator=()

G4PenelopeBremsstrahlungModel & G4PenelopeBremsstrahlungModel::operator= ( const G4PenelopeBremsstrahlungModel right)
delete

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 322 of file G4PenelopeBremsstrahlungModel.cc.

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

103{fVerboseLevel = 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: