Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
GB05BOptrSplitAndKillByCrossSection Class Reference

#include <Doxymodules_biasing.h>

Inheritance diagram for GB05BOptrSplitAndKillByCrossSection:
G4VBiasingOperator

Public Member Functions

 GB05BOptrSplitAndKillByCrossSection (G4String particleToBias, G4String name="SplitAndKillByXS")
 
virtual ~GB05BOptrSplitAndKillByCrossSection ()
 
virtual void StartRun ()
 
void AddProcessToEquipoise (G4String processName)
 

Private Member Functions

virtual G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *, const G4BiasingProcessInterface *) final
 
virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *, const G4BiasingProcessInterface *) final
 
virtual G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *, const G4BiasingProcessInterface *) final
 

Private Attributes

GB05BOptnSplitAndKillByCrossSectionfSplitAndKillByCrossSection
 
const G4ParticleDefinitionfParticleToBias
 
std::vector< G4StringfProcessesToEquipoise
 
G4bool fSetup
 
std::vector< const G4VProcess * > fProcesses
 

Detailed Description

Definition at line 124 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ GB05BOptrSplitAndKillByCrossSection()

GB05BOptrSplitAndKillByCrossSection::GB05BOptrSplitAndKillByCrossSection ( G4String  particleToBias,
G4String  name = "SplitAndKillByXS" 
)

Definition at line 40 of file GB05BOptrSplitAndKillByCrossSection.cc.

43 : G4VBiasingOperator(name),
44 fSetup(true)
45{
46 fParticleToBias = G4ParticleTable::GetParticleTable()->FindParticle(particleName);
47
48 if ( fParticleToBias == 0 )
49 {
50 G4ExceptionDescription ed;
51 ed << "Particle `" << particleName << "' not found !" << G4endl;
52 G4Exception("GB05BOptrSplitAndKillByCrossSection(...)",
53 "exGB05.01",
54 JustWarning,
55 ed);
56 }
57
59 new GB05BOptnSplitAndKillByCrossSection("splitterFor_"+particleName);
60
61}
GB05BOptnSplitAndKillByCrossSection * fSplitAndKillByCrossSection

◆ ~GB05BOptrSplitAndKillByCrossSection()

GB05BOptrSplitAndKillByCrossSection::~GB05BOptrSplitAndKillByCrossSection ( )
virtual

Definition at line 65 of file GB05BOptrSplitAndKillByCrossSection.cc.

66{
68}

Member Function Documentation

◆ StartRun()

void GB05BOptrSplitAndKillByCrossSection::StartRun ( )
virtual

Definition at line 72 of file GB05BOptrSplitAndKillByCrossSection.cc.

73{
74 // ---------------
75 // -- Setup stage:
76 // ---------------
77 // -- Start by collecting the pointer of the physics processes
78 // -- considered for the splitting by cross-sections. Doing so,
79 // -- this also verifies that these physics processes are each
80 // -- under control of a G4BiasingProcessInterface wrapper.
81 if ( fSetup )
82 {
83 const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
84 const G4BiasingProcessSharedData* sharedData =
85 G4BiasingProcessInterface::GetSharedData( processManager );
86 if ( sharedData )
87 {
88 for ( size_t i = 0 ; i < fProcessesToEquipoise.size() ; i++ )
89 {
90 G4bool processFound(false);
91 for ( size_t j = 0 ;
92 j < (sharedData->GetPhysicsBiasingProcessInterfaces()).size();
93 j++ )
94 {
95 const G4BiasingProcessInterface* wrapperProcess =
96 (sharedData->GetPhysicsBiasingProcessInterfaces())[j];
97 if ( fProcessesToEquipoise[i] ==
98 wrapperProcess->GetWrappedProcess()->GetProcessName() )
99 {
100 fProcesses.push_back( wrapperProcess->GetWrappedProcess() );
101 processFound = true;
102 break;
103 }
104 }
105 if ( !processFound )
106 {
107 G4String particleName = "(unknown)";
108 if ( fParticleToBias != nullptr )
109 {
110 particleName = fParticleToBias->GetParticleName();
111 }
112 G4ExceptionDescription ed;
113 ed << "Process `" << fProcessesToEquipoise[i]
114 << "' not found for particle `" << particleName << "'"
115 << G4endl;
116 G4Exception("GB05BOptrSplitAndKillByCrossSection::StartRun(...)",
117 "exGB05.02",
118 JustWarning,
119 ed);
120 }
121 }
122 }
123 fSetup = false;
124 }
125
126 if ( fProcessesToEquipoise.size() == 0 || fProcesses.size() == 0 )
127 {
128 G4ExceptionDescription ed;
129 ed << "No processes to counterbalance for defined or found ! "
130 << "Biasing will do nothing."
131 << G4endl;
132 G4Exception("GB05BOptrSplitAndKillByCrossSection::StartRun(...)",
133 "exGB05.03",
134 JustWarning,
135 ed);
136 }
137
138}

◆ ProposeOccurenceBiasingOperation()

virtual G4VBiasingOperation * GB05BOptrSplitAndKillByCrossSection::ProposeOccurenceBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlinefinalprivatevirtual

Definition at line 75 of file GB05BOptrSplitAndKillByCrossSection.hh.

77 { return 0; }

◆ ProposeFinalStateBiasingOperation()

virtual G4VBiasingOperation * GB05BOptrSplitAndKillByCrossSection::ProposeFinalStateBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlinefinalprivatevirtual

Definition at line 81 of file GB05BOptrSplitAndKillByCrossSection.hh.

83 { return 0; }

◆ ProposeNonPhysicsBiasingOperation()

G4VBiasingOperation * GB05BOptrSplitAndKillByCrossSection::ProposeNonPhysicsBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface  
)
finalprivatevirtual

Definition at line 143 of file GB05BOptrSplitAndKillByCrossSection.cc.

146{
147
148 // -----------------------------------------------------
149 // -- Check if current particle type is the one to bias:
150 // -----------------------------------------------------
151 if ( track->GetDefinition() != fParticleToBias ) return 0;
152
153 // --------------------------------------------------------------------
154 // -- Compute the total cross-section for the physics processes
155 // -- considered.
156 // -- These physics processes have been updated to the current track
157 // -- state by their related wrapper G4BiasingProcessInterface objects,
158 // -- the cross-sections/mean free pathes are hence usable.
159 // --------------------------------------------------------------------
160 G4double totalCrossSection(0.0);
161 for ( size_t i = 0 ; i < fProcesses.size() ; i++ )
162 {
163 G4double interactionLength = fProcesses[i]->GetCurrentInteractionLength();
164 if ( interactionLength < DBL_MAX/10. )
165 totalCrossSection += 1./interactionLength;
166 }
167 if ( totalCrossSection < DBL_MIN ) return 0;
168
169 G4double totalInteractionLength = 1./totalCrossSection;
170
171 // ---------------------------------------------------------------------
172 // -- Passes the updated "absorption" cross-section (interaction length)
173 // -- to the biasing operation, and returns this operation:
174 // ---------------------------------------------------------------------
175 fSplitAndKillByCrossSection->SetInteractionLength( totalInteractionLength );
176
178
179}

◆ AddProcessToEquipoise()

void GB05BOptrSplitAndKillByCrossSection::AddProcessToEquipoise ( G4String  processName)

Definition at line 181 of file GB05BOptrSplitAndKillByCrossSection.cc.

183{
184 fProcessesToEquipoise.push_back( processName );
185}

Member Data Documentation

◆ fSplitAndKillByCrossSection

GB05BOptnSplitAndKillByCrossSection* GB05BOptrSplitAndKillByCrossSection::fSplitAndKillByCrossSection
private

Definition at line 102 of file GB05BOptrSplitAndKillByCrossSection.hh.

◆ fParticleToBias

const G4ParticleDefinition* GB05BOptrSplitAndKillByCrossSection::fParticleToBias
private

Definition at line 103 of file GB05BOptrSplitAndKillByCrossSection.hh.

◆ fProcessesToEquipoise

std::vector< G4String > GB05BOptrSplitAndKillByCrossSection::fProcessesToEquipoise
private

Definition at line 104 of file GB05BOptrSplitAndKillByCrossSection.hh.

◆ fSetup

G4bool GB05BOptrSplitAndKillByCrossSection::fSetup
private

Definition at line 105 of file GB05BOptrSplitAndKillByCrossSection.hh.

◆ fProcesses

std::vector< const G4VProcess* > GB05BOptrSplitAndKillByCrossSection::fProcesses
private

Definition at line 106 of file GB05BOptrSplitAndKillByCrossSection.hh.


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

Applications | User Support | Publications | Collaboration