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

#include <G4ChannelingOptrChangeCrossSection.hh>

+ Inheritance diagram for G4ChannelingOptrChangeCrossSection:

Public Member Functions

 G4ChannelingOptrChangeCrossSection (G4String particleToBias, G4String name="ChannelingChangeXS")
 
virtual ~G4ChannelingOptrChangeCrossSection ()
 
virtual void StartRun ()
 
- Public Member Functions inherited from G4VBiasingOperator
 G4VBiasingOperator (G4String name)
 
virtual ~G4VBiasingOperator ()
 
virtual void Configure ()
 
virtual void ConfigureForWorker ()
 
virtual void StartRun ()
 
virtual void StartTracking (const G4Track *)
 
virtual void EndTracking ()
 
const G4String GetName () const
 
void AttachTo (const G4LogicalVolume *)
 
G4BiasingAppliedCase GetPreviousBiasingAppliedCase () const
 
G4VBiasingOperationGetProposedOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
void ExitingBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
const G4VBiasingOperationGetPreviousNonPhysicsAppliedOperation ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VBiasingOperator
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators ()
 
static G4VBiasingOperatorGetBiasingOperator (const G4LogicalVolume *)
 
- Protected Member Functions inherited from G4VBiasingOperator
virtual G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void ExitBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 

Detailed Description

Definition at line 60 of file G4ChannelingOptrChangeCrossSection.hh.

Constructor & Destructor Documentation

◆ G4ChannelingOptrChangeCrossSection()

G4ChannelingOptrChangeCrossSection::G4ChannelingOptrChangeCrossSection ( G4String  particleToBias,
G4String  name = "ChannelingChangeXS" 
)

Definition at line 43 of file G4ChannelingOptrChangeCrossSection.cc.

46fChannelingID(-1),
47fSetup(true){
48 fParticleToBias = G4ParticleTable::GetParticleTable()->FindParticle(particleName);
49
50 if ( fParticleToBias == 0 )
51 {
53 ed << "Particle `" << particleName << "' not found !" << G4endl;
54 G4Exception("G4ChannelingOptrChangeCrossSection(...)",
55 "G4Channeling",
57 ed);
58 }
59
60 fProcessToDensity["channeling"] = fDensityRatioNone;
61}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4endl
Definition: G4ios.hh:57
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()

◆ ~G4ChannelingOptrChangeCrossSection()

G4ChannelingOptrChangeCrossSection::~G4ChannelingOptrChangeCrossSection ( )
virtual

Definition at line 65 of file G4ChannelingOptrChangeCrossSection.cc.

65 {
66 for ( std::map< const G4BiasingProcessInterface*, G4BOptnChangeCrossSection* >::iterator
67 it = fChangeCrossSectionOperations.begin() ;
68 it != fChangeCrossSectionOperations.end() ;
69 it++ ) delete (*it).second;
70}

Member Function Documentation

◆ StartRun()

void G4ChannelingOptrChangeCrossSection::StartRun ( )
virtual

Reimplemented from G4VBiasingOperator.

Definition at line 74 of file G4ChannelingOptrChangeCrossSection.cc.

74 {
75 if ( fSetup ){
76 const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
77 const G4BiasingProcessSharedData* sharedData =
79 if ( sharedData ){
80 for ( size_t i = 0 ; i < (sharedData->GetPhysicsBiasingProcessInterfaces()).size(); i++ ){
81 const G4BiasingProcessInterface* wrapperProcess =
82 (sharedData->GetPhysicsBiasingProcessInterfaces())[i];
83 G4String processName = wrapperProcess->GetWrappedProcess()->GetProcessName();
84 G4String operationName = "channelingChangeXS-" + processName;
85 fChangeCrossSectionOperations[wrapperProcess] =
86 new G4BOptnChangeCrossSection(operationName);
87
88 G4ProcessType type = wrapperProcess->GetWrappedProcess()->GetProcessType();
89 G4int subType = wrapperProcess->GetWrappedProcess()->GetProcessSubType();
90
91 switch (type) {
92 case fNotDefined:
93 fProcessToDensity[processName] = fDensityRatioNotDefined;
94 break;
95 case fTransportation:
96 fProcessToDensity[processName] = fDensityRatioNone;
97 break;
99 if(subType == fCoulombScattering ||
100 subType == fMultipleScattering){
101 fProcessToDensity[processName] = fDensityRatioNuD;
102 }
103 if(subType == fIonisation ||
104 subType == fPairProdByCharged ||
105 subType == fAnnihilation ||
106 subType == fAnnihilationToMuMu ||
107 subType == fAnnihilationToHadrons){
108 fProcessToDensity[processName] = fDensityRatioElD;
109 }
110 if(subType == fBremsstrahlung ||
111 subType == fNuclearStopping){
112 fProcessToDensity[processName] = fDensityRatioNuDElD;
113 }
114
115 if(subType == fCerenkov ||
116 subType == fScintillation ||
117 subType == fSynchrotronRadiation ||
118 subType == fTransitionRadiation){
119 fProcessToDensity[processName] = fDensityRatioNone;
120 }
121 if(subType == fRayleigh ||
122 subType == fPhotoElectricEffect ||
123 subType == fComptonScattering ||
124 subType == fGammaConversion ||
125 subType == fGammaConversionToMuMu){
126 fProcessToDensity[processName] = fDensityRatioNone;
127 }
128 break;
129 case fOptical:
130 fProcessToDensity[processName] = fDensityRatioNone;
131 break;
132 case fHadronic:
133 fProcessToDensity[processName] = fDensityRatioNuD;
134 break;
136 fProcessToDensity[processName] = fDensityRatioNuD;
137 break;
138 case fGeneral:
139 fProcessToDensity[processName] = fDensityRatioNone;
140 break;
141 case fDecay:
142 fProcessToDensity[processName] = fDensityRatioNone;
143 break;
145 fProcessToDensity[processName] = fDensityRatioNone;
146 break;
147 case fUserDefined:
148 fProcessToDensity[processName] = fDensityRatioNone;
149 break;
150 case fParallel:
151 fProcessToDensity[processName] = fDensityRatioNone;
152 break;
153 case fPhonon:
154 fProcessToDensity[processName] = fDensityRatioNone;
155 break;
156 case fUCN:
157 fProcessToDensity[processName] = fDensityRatioNone;
158 break;
159 default:
160 fProcessToDensity[processName] = fDensityRatioNone;
161 break;
162 }
163 }
164 }
165 fSetup = false;
166 }
167}
@ fGammaConversionToMuMu
@ fAnnihilationToHadrons
@ fBremsstrahlung
@ fCoulombScattering
@ fGammaConversion
@ fRayleigh
@ fIonisation
@ fPairProdByCharged
@ fSynchrotronRadiation
@ fCerenkov
@ fAnnihilationToMuMu
@ fScintillation
@ fNuclearStopping
@ fComptonScattering
@ fTransitionRadiation
@ fAnnihilation
@ fMultipleScattering
@ fPhotoElectricEffect
G4ProcessType
@ fOptical
@ fPhonon
@ fParameterisation
@ fParallel
@ fUCN
@ fGeneral
@ fDecay
@ fElectromagnetic
@ fHadronic
@ fUserDefined
@ fTransportation
@ fPhotolepton_hadron
@ fNotDefined
int G4int
Definition: G4Types.hh:85
const G4BiasingProcessSharedData * GetSharedData() const
G4VProcess * GetWrappedProcess() const
const std::vector< const G4BiasingProcessInterface * > & GetPhysicsBiasingProcessInterfaces() const
G4ProcessManager * GetProcessManager() const
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:388
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

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