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

#include <G4CoherentPairProductionPhysics.hh>

+ Inheritance diagram for G4CoherentPairProductionPhysics:

Public Member Functions

 G4CoherentPairProductionPhysics (const G4String &name="Coherent Pair Production Physics")
 
 ~G4CoherentPairProductionPhysics ()=default
 
void ConstructParticle () override
 
void ConstructProcess () override
 
void ActivateIncoherentScattering ()
 
void SetNameChannelingModel (const G4String &nameChannelingModel)
 set functions
 
void SetNameG4Region (const G4String &nameG4Region)
 set name of G4Region to where the G4ChannelingFastSimModel is active
 
void SetLowEnergyLimit (G4double energy)
 
void SetHighAngleLimit (G4double angle)
 
void SetPPKineticEnergyCut (G4double kineticEnergyCut)
 
void SetSamplingPairsNumber (G4int nPairs)
 
void SetChargeParticleAngleFactor (G4double chargeParticleAngleFactor)
 
void SetNTrajectorySteps (G4int nTrajectorySteps)
 set number of trajectory steps of a single particle (e- or e+)
 
- Public Member Functions inherited from G4VPhysicsConstructor
 G4VPhysicsConstructor (const G4String &="")
 
 G4VPhysicsConstructor (const G4String &name, G4int physics_type)
 
virtual ~G4VPhysicsConstructor ()
 
void SetPhysicsName (const G4String &="")
 
const G4StringGetPhysicsName () const
 
void SetPhysicsType (G4int)
 
G4int GetPhysicsType () const
 
G4int GetInstanceID () const
 
virtual void TerminateWorker ()
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VPhysicsConstructor
static const G4VPCManagerGetSubInstanceManager ()
 
- Protected Types inherited from G4VPhysicsConstructor
using PhysicsBuilder_V = G4VPCData::PhysicsBuilders_V
 
- Protected Member Functions inherited from G4VPhysicsConstructor
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 
PhysicsBuilder_V GetBuilders () const
 
void AddBuilder (G4PhysicsBuilderInterface *bld)
 
- Protected Attributes inherited from G4VPhysicsConstructor
G4int verboseLevel = 0
 
G4String namePhysics = ""
 
G4int typePhysics = 0
 
G4ParticleTabletheParticleTable = nullptr
 
G4int g4vpcInstanceID = 0
 
- Static Protected Attributes inherited from G4VPhysicsConstructor
static G4RUN_DLL G4VPCManager subInstanceManager
 

Detailed Description

Definition at line 36 of file G4CoherentPairProductionPhysics.hh.

Constructor & Destructor Documentation

◆ G4CoherentPairProductionPhysics()

G4CoherentPairProductionPhysics::G4CoherentPairProductionPhysics ( const G4String & name = "Coherent Pair Production Physics")

Definition at line 38 of file G4CoherentPairProductionPhysics.cc.

38 :
40{}
G4VPhysicsConstructor(const G4String &="")

◆ ~G4CoherentPairProductionPhysics()

G4CoherentPairProductionPhysics::~G4CoherentPairProductionPhysics ( )
default

Member Function Documentation

◆ ActivateIncoherentScattering()

void G4CoherentPairProductionPhysics::ActivateIncoherentScattering ( )
inline

activate incoherent scattering (standard gamma conversion should be switched off in physics list)

Definition at line 49 of file G4CoherentPairProductionPhysics.hh.

49{fIncoherentScattering = true;}

◆ ConstructParticle()

void G4CoherentPairProductionPhysics::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 44 of file G4CoherentPairProductionPhysics.cc.

45{
46 G4Gamma::GammaDefinition(); // Define the gamma particle
47}
static G4Gamma * GammaDefinition()
Definition G4Gamma.cc:76

◆ ConstructProcess()

void G4CoherentPairProductionPhysics::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 51 of file G4CoherentPairProductionPhysics.cc.

52{
53 //create the gamma process
54 G4CoherentPairProduction* gammaProcess = new G4CoherentPairProduction();
55
56 if(verboseLevel > 0) {
57 G4cout << "G4CoherentPairProductionPhysics::ConstructProcess" << G4endl;
58 }
59
60 G4RegionStore* regionStore = G4RegionStore::GetInstance();
61 G4Region* RegionCh = regionStore->GetRegion(fNameRegion);
62
63 //if the region is not found
64 if(RegionCh==0)
65 {
66 G4Exception("GetRegion",// Origin of the exception
67 "001", // Unique error code
68 FatalException, // Terminate the program
69 "Region is not found! The program will terminate.");
70 }
71 else
72 {
73 //get channeling model
74 G4bool someflag=false;
75 G4ChannelingFastSimModel* Channeling =
76 static_cast<G4ChannelingFastSimModel*>
77 (RegionCh->GetFastSimulationManager()->
78 GetFastSimulationModel(fNameChannelingModel,0,someflag));
79
80 //if channeling model is not found
81 if(Channeling==0)
82 {
83 G4Exception("GetFastSimulationModel",// Origin of the exception
84 "001", // Unique error code
85 FatalException, // Terminate the program
86 "Input channeling model is not found! The program will terminate."
87 );
88 }
89 else
90 {
91 gammaProcess->Input(Channeling->GetCrystalData());
92 }
93 }
94
95 //set functions
96 if(fIncoherentScattering){gammaProcess->ActivateIncoherentScattering();}
97 gammaProcess->SetLowEnergyLimit(fLowEnergyLimit);
98 gammaProcess->SetHighAngleLimit(fHighAngleLimit);
99 gammaProcess->SetPPKineticEnergyCut(fPPKineticEnergyCut);
100 gammaProcess->SetSamplingPairsNumber(fNMCPairs);
101 gammaProcess->SetChargeParticleAngleFactor(fChargeParticleAngleFactor);
102 gammaProcess->SetNTrajectorySteps(fNTrajectorySteps);
103
104 //set the name of G4Region in which the model is applicable
105 gammaProcess->SetG4RegionName(fNameRegion);
106
107 //get process manager for gamma
108 G4ProcessManager* pManager = G4Gamma::Gamma()->GetProcessManager();
109
110 //register the G4CoherentPairProduction process
111 pManager->AddDiscreteProcess(gammaProcess);
112}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ChannelingFastSimCrystalData * GetCrystalData()
void SetPPKineticEnergyCut(G4double kineticEnergyCut)
void SetLowEnergyLimit(G4double energy)
set cuts
void SetG4RegionName(const G4String &nameG4Region)
set the name of G4Region in which the model is applicable
void Input(const G4Material *crystal, const G4String &lattice)
special functions
void SetNTrajectorySteps(G4int nTrajectorySteps)
set number of trajectory steps of a single particle (e- or e+)
void SetChargeParticleAngleFactor(G4double chargeParticleAngleFactor)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4ProcessManager * GetProcessManager() const
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
G4FastSimulationManager * GetFastSimulationManager() const
Definition G4Region.cc:140

◆ SetChargeParticleAngleFactor()

void G4CoherentPairProductionPhysics::SetChargeParticleAngleFactor ( G4double chargeParticleAngleFactor)
inline

set the number of particle angles 1/gamma in pair production defining the width of the angular distribution of pair sampling in the Baier-Katkov Integral

Definition at line 72 of file G4CoherentPairProductionPhysics.hh.

73 {fChargeParticleAngleFactor = chargeParticleAngleFactor;}

◆ SetHighAngleLimit()

void G4CoherentPairProductionPhysics::SetHighAngleLimit ( G4double angle)
inline

Definition at line 62 of file G4CoherentPairProductionPhysics.hh.

62{fHighAngleLimit=angle;}

◆ SetLowEnergyLimit()

void G4CoherentPairProductionPhysics::SetLowEnergyLimit ( G4double energy)
inline

Definition at line 61 of file G4CoherentPairProductionPhysics.hh.

61{fLowEnergyLimit=energy;}
G4double energy(const ThreeVector &p, const G4double m)

◆ SetNameChannelingModel()

void G4CoherentPairProductionPhysics::SetNameChannelingModel ( const G4String & nameChannelingModel)
inline

set functions

set name of G4ChannelingFastSimModel from which an auto input should be performed

Definition at line 54 of file G4CoherentPairProductionPhysics.hh.

55 {fNameChannelingModel=nameChannelingModel;}

◆ SetNameG4Region()

void G4CoherentPairProductionPhysics::SetNameG4Region ( const G4String & nameG4Region)
inline

set name of G4Region to where the G4ChannelingFastSimModel is active

Definition at line 58 of file G4CoherentPairProductionPhysics.hh.

59 {fNameRegion=nameG4Region;}

◆ SetNTrajectorySteps()

void G4CoherentPairProductionPhysics::SetNTrajectorySteps ( G4int nTrajectorySteps)
inline

set number of trajectory steps of a single particle (e- or e+)

Definition at line 76 of file G4CoherentPairProductionPhysics.hh.

77 {fNTrajectorySteps = nTrajectorySteps;}

◆ SetPPKineticEnergyCut()

void G4CoherentPairProductionPhysics::SetPPKineticEnergyCut ( G4double kineticEnergyCut)
inline

Definition at line 63 of file G4CoherentPairProductionPhysics.hh.

63{fPPKineticEnergyCut=kineticEnergyCut;}

◆ SetSamplingPairsNumber()

void G4CoherentPairProductionPhysics::SetSamplingPairsNumber ( G4int nPairs)
inline

set the number of pairs in sampling of Baier-Katkov Integral (MC integration by e+- energy and angles <=> e+- momentum)

Definition at line 67 of file G4CoherentPairProductionPhysics.hh.

67{fNMCPairs = nPairs;}

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