Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PenelopeBremsstrahlungModel Class Reference

#include <G4PenelopeBremsstrahlungModel.hh>

Inheritance diagram for G4PenelopeBremsstrahlungModel:
Collaboration 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 InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

Protected Attributes

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

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 *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 63 of file G4PenelopeBremsstrahlungModel.hh.

Constructor & Destructor Documentation

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

Definition at line 64 of file G4PenelopeBremsstrahlungModel.cc.

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

Here is the call graph for this function:

G4PenelopeBremsstrahlungModel::~G4PenelopeBremsstrahlungModel ( )
virtual

Definition at line 99 of file G4PenelopeBremsstrahlungModel.cc.

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

Here is the call graph for this function:

Member Function Documentation

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

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

Reimplemented from G4VEmModel.

Definition at line 298 of file G4PenelopeBremsstrahlungModel.cc.

302 {
303  if (verboseLevel > 3)
304  G4cout << "Calling ComputeDEDX() of G4PenelopeBremsstrahlungModel" << G4endl;
305 
306  G4PenelopeCrossSection* theXS = GetCrossSectionTableForCouple(theParticle,material,
307  cutEnergy);
308  G4double sPowerPerMolecule = 0.0;
309  if (theXS)
310  sPowerPerMolecule = theXS->GetSoftStoppingPower(kineticEnergy);
311 
312 
313  G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
314  G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
315 
316  G4double moleculeDensity = 0.;
317  if (atPerMol)
318  moleculeDensity = atomDensity/atPerMol;
319 
320  G4double sPowerPerVolume = sPowerPerMolecule*moleculeDensity;
321 
322  if (verboseLevel > 2)
323  {
324  G4cout << "G4PenelopeBremsstrahlungModel " << G4endl;
325  G4cout << "Stopping power < " << cutEnergy/keV << " keV at " <<
326  kineticEnergy/keV << " keV = " <<
327  sPowerPerVolume/(keV/mm) << " keV/mm" << G4endl;
328  }
329  return sPowerPerVolume;
330 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetSoftStoppingPower(G4double energy) const
Returns the total stopping power due to soft collisions.
G4GLOB_DLL std::ostream G4cout
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:209
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetAtomsPerMolecule(const G4Material *)
Returns the total number of atoms per molecule.
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 234 of file G4PenelopeBremsstrahlungModel.cc.

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

Here is the call graph for this function:

G4int G4PenelopeBremsstrahlungModel::GetVerbosityLevel ( )
inline

Definition at line 109 of file G4PenelopeBremsstrahlungModel.hh.

109 {return verboseLevel;};
void G4PenelopeBremsstrahlungModel::Initialise ( const G4ParticleDefinition part,
const G4DataVector theCuts 
)
virtual

Implements G4VEmModel.

Definition at line 114 of file G4PenelopeBremsstrahlungModel.cc.

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

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 180 of file G4PenelopeBremsstrahlungModel.cc.

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

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 473 of file G4PenelopeBremsstrahlungModel.cc.

475 {
476  return fIntrinsicLowEnergyLimit;
477 }
void G4PenelopeBremsstrahlungModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicParticle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 334 of file G4PenelopeBremsstrahlungModel.cc.

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

Here is the call graph for this function:

void G4PenelopeBremsstrahlungModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 108 of file G4PenelopeBremsstrahlungModel.hh.

108 {verboseLevel = lev;};

Member Data Documentation

const G4ParticleDefinition* G4PenelopeBremsstrahlungModel::fParticle
protected

Definition at line 113 of file G4PenelopeBremsstrahlungModel.hh.

G4ParticleChangeForLoss* G4PenelopeBremsstrahlungModel::fParticleChange
protected

Definition at line 109 of file G4PenelopeBremsstrahlungModel.hh.


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