Geant4 11.1.1
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 44 of file G4ChannelingOptrChangeCrossSection.cc.

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

◆ ~G4ChannelingOptrChangeCrossSection()

G4ChannelingOptrChangeCrossSection::~G4ChannelingOptrChangeCrossSection ( )
virtual

Definition at line 66 of file G4ChannelingOptrChangeCrossSection.cc.

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

Member Function Documentation

◆ StartRun()

void G4ChannelingOptrChangeCrossSection::StartRun ( )
virtual

Reimplemented from G4VBiasingOperator.

Definition at line 75 of file G4ChannelingOptrChangeCrossSection.cc.

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

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