Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointBremsstrahlungModel.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26////////////////////////////////////////////////////////////////////////////////
27// Class: G4AdjointBremsstrahlungModel
28// Author: L. Desorgher
29// Organisation: SpaceIT GmbH
30////////////////////////////////////////////////////////////////////////////////
31//
32// Adjoint Model for e- Bremsstrahlung.Adapted from G4eBremsstrahlungModel
33// Use of a simple biased differential cross section (C(Z)/Egamma) allowing a
34// rapid computation of adjoint CS and rapid sampling of adjoint secondaries.
35// In this way cross section matrices are not used anymore, avoiding a long
36// computation of adjoint brem cross section matrices for each material
37// at initialisation. This mode can be switched on/off by selecting
38// SetUseMatrix(false)/ SetUseMatrix(true) in the constructor.
39//
40//-------------------------------------------------------------
41
42#ifndef G4AdjointBremsstrahlungModel_h
43#define G4AdjointBremsstrahlungModel_h 1
44
45#include "globals.hh"
46#include "G4VEmAdjointModel.hh"
47
51
53{
54 public:
56
58
60
61 void SampleSecondaries(const G4Track& aTrack, G4bool isScatProjToProj,
62 G4ParticleChange* fParticleChange) override;
63
64 void RapidSampleSecondaries(const G4Track& aTrack, G4bool isScatProjToProj,
65 G4ParticleChange* fParticleChange);
66
68 const G4Material* aMaterial,
69 G4double kinEnergyProj, // kinetic energy of the primary particle before
70 // the interaction
71 G4double kinEnergyProd // kinetic energy of the secondary particle
72 ) override;
73
75 G4double primEnergy,
76 G4bool isScatProjToProj) override;
77
80 const G4AdjointBremsstrahlungModel& right) = delete;
81
82 private:
83 void Initialize();
84
85 G4EmModelManager* fEmModelManagerForFwdModels;
86 G4AdjointCSManager* fCSManager;
87 G4ParticleDefinition* fElectron;
89
90 G4double fLastCZ = 0.;
91
92 G4bool fIsDirectModelInitialised = false;
93};
94
95#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
void SampleSecondaries(const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange) override
G4double DiffCrossSectionPerVolumePrimToSecond(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd) override
G4AdjointBremsstrahlungModel & operator=(const G4AdjointBremsstrahlungModel &right)=delete
G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool isScatProjToProj) override
void RapidSampleSecondaries(const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange)
G4AdjointBremsstrahlungModel(G4AdjointBremsstrahlungModel &)=delete