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

#include <Doxymodules_biasing.h>

Inheritance diagram for GB01BOptrChangeCrossSection:
G4VBiasingOperator

Public Member Functions

 GB01BOptrChangeCrossSection (G4String particleToBias, G4String name="ChangeXS")
 
virtual ~GB01BOptrChangeCrossSection ()
 
virtual void StartRun ()
 

Private Member Functions

virtual G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *, const G4BiasingProcessInterface *)
 
virtual G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *, const G4BiasingProcessInterface *)
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 

Private Attributes

std::map< const G4BiasingProcessInterface *, G4BOptnChangeCrossSection * > fChangeCrossSectionOperations
 
G4bool fSetup
 
const G4ParticleDefinitionfParticleToBias
 

Detailed Description

Definition at line 67 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ GB01BOptrChangeCrossSection()

GB01BOptrChangeCrossSection::GB01BOptrChangeCrossSection ( G4String  particleToBias,
G4String  name = "ChangeXS" 
)

Definition at line 43 of file GB01BOptrChangeCrossSection.cc.

45 : G4VBiasingOperator(name),
46 fSetup(true)
47{
48 fParticleToBias = G4ParticleTable::GetParticleTable()->FindParticle(particleName);
49
50 if ( fParticleToBias == 0 )
51 {
52 G4ExceptionDescription ed;
53 ed << "Particle `" << particleName << "' not found !" << G4endl;
54 G4Exception("GB01BOptrChangeCrossSection(...)",
55 "exGB01.01",
56 JustWarning,
57 ed);
58 }
59}
const G4ParticleDefinition * fParticleToBias

◆ ~GB01BOptrChangeCrossSection()

GB01BOptrChangeCrossSection::~GB01BOptrChangeCrossSection ( )
virtual

Definition at line 63 of file GB01BOptrChangeCrossSection.cc.

64{
65 for ( std::map< const G4BiasingProcessInterface*, G4BOptnChangeCrossSection* >::iterator
68 it++ ) delete (*it).second;
69}
std::map< const G4BiasingProcessInterface *, G4BOptnChangeCrossSection * > fChangeCrossSectionOperations

Member Function Documentation

◆ StartRun()

void GB01BOptrChangeCrossSection::StartRun ( )
virtual

Definition at line 73 of file GB01BOptrChangeCrossSection.cc.

74{
75 // --------------
76 // -- Setup stage:
77 // ---------------
78 // -- Start by collecting processes under biasing, create needed biasing
79 // -- operations and associate these operations to the processes:
80 if ( fSetup )
81 {
82 const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
83 const G4BiasingProcessSharedData* sharedData =
84 G4BiasingProcessInterface::GetSharedData( processManager );
85 if ( sharedData ) // -- sharedData tested, as is can happen a user attaches an operator to a
86 // -- volume but without defined BiasingProcessInterface processes.
87 {
88 for ( size_t i = 0 ; i < (sharedData->GetPhysicsBiasingProcessInterfaces()).size(); i++ )
89 {
90 const G4BiasingProcessInterface* wrapperProcess =
91 (sharedData->GetPhysicsBiasingProcessInterfaces())[i];
92 G4String operationName = "XSchange-" +
93 wrapperProcess->GetWrappedProcess()->GetProcessName();
94 fChangeCrossSectionOperations[wrapperProcess] =
95 new G4BOptnChangeCrossSection(operationName);
96 }
97 }
98 fSetup = false;
99 }
100}

◆ ProposeOccurenceBiasingOperation()

G4VBiasingOperation * GB01BOptrChangeCrossSection::ProposeOccurenceBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
privatevirtual

Definition at line 105 of file GB01BOptrChangeCrossSection.cc.

108{
109
110 // -----------------------------------------------------
111 // -- Check if current particle type is the one to bias:
112 // -----------------------------------------------------
113 if ( track->GetDefinition() != fParticleToBias ) return 0;
114
115
116 // ---------------------------------------------------------------------
117 // -- select and setup the biasing operation for current callingProcess:
118 // ---------------------------------------------------------------------
119 // -- Check if the analog cross-section well defined : for example, the conversion
120 // -- process for a gamma below e+e- creation threshold has an DBL_MAX interaction
121 // -- length. Nothing is done in this case (ie, let analog process to deal with the case)
122 G4double analogInteractionLength =
123 callingProcess->GetWrappedProcess()->GetCurrentInteractionLength();
124 if ( analogInteractionLength > DBL_MAX/10. ) return 0;
125
126 // -- Analog cross-section is well-defined:
127 G4double analogXS = 1./analogInteractionLength;
128
129 // -- Choose a constant cross-section bias. But at this level, this factor can be made
130 // -- direction dependent, like in the exponential transform MCNP case, or it
131 // -- can be chosen differently, depending on the process, etc.
132 G4double XStransformation = 2.0 ;
133
134 // -- fetch the operation associated to this callingProcess:
135 G4BOptnChangeCrossSection* operation = fChangeCrossSectionOperations[callingProcess];
136 // -- get the operation that was proposed to the process in the previous step:
137 G4VBiasingOperation* previousOperation = callingProcess->GetPreviousOccurenceBiasingOperation();
138
139 // -- now setup the operation to be returned to the process: this
140 // -- consists in setting the biased cross-section, and in asking
141 // -- the operation to sample its exponential interaction law.
142 // -- To do this, to first order, the two lines:
143 // operation->SetBiasedCrossSection( XStransformation * analogXS );
144 // operation->Sample();
145 // -- are correct and sufficient.
146 // -- But, to avoid having to shoot a random number each time, we sample
147 // -- only on the first time the operation is proposed, or if the interaction
148 // -- occured. If the interaction did not occur for the process in the previous,
149 // -- we update the number of interaction length instead of resampling.
150 if ( previousOperation == 0 )
151 {
152 operation->SetBiasedCrossSection( XStransformation * analogXS );
153 operation->Sample();
154 }
155 else
156 {
157 if ( previousOperation != operation )
158 {
159 // -- should not happen !
160 G4ExceptionDescription ed;
161 ed << " Logic problem in operation handling !" << G4endl;
162 G4Exception("GB01BOptrChangeCrossSection::ProposeOccurenceBiasingOperation(...)",
163 "exGB01.02",
164 JustWarning,
165 ed);
166 return 0;
167 }
168 if ( operation->GetInteractionOccured() )
169 {
170 operation->SetBiasedCrossSection( XStransformation * analogXS );
171 operation->Sample();
172 }
173 else
174 {
175 // -- update the 'interaction length' and underneath 'number of interaction lengths'
176 // -- for past step (this takes into accout the previous step cross-section value)
177 operation->UpdateForStep( callingProcess->GetPreviousStepSize() );
178 // -- update the cross-section value:
179 operation->SetBiasedCrossSection( XStransformation * analogXS );
180 // -- forces recomputation of the 'interaction length' taking into account above
181 // -- new cross-section value [tricky : to be improved]
182 operation->UpdateForStep( 0.0 );
183 }
184 }
185
186 return operation;
187
188}

◆ ProposeFinalStateBiasingOperation()

virtual G4VBiasingOperation * GB01BOptrChangeCrossSection::ProposeFinalStateBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlineprivatevirtual

Definition at line 75 of file GB01BOptrChangeCrossSection.hh.

76 {return 0;}

◆ ProposeNonPhysicsBiasingOperation()

virtual G4VBiasingOperation * GB01BOptrChangeCrossSection::ProposeNonPhysicsBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlineprivatevirtual

Definition at line 78 of file GB01BOptrChangeCrossSection.hh.

79 {return 0;}

◆ OperationApplied()

void GB01BOptrChangeCrossSection::OperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation occurenceOperationApplied,
G4double  weightForOccurenceInteraction,
G4VBiasingOperation finalStateOperationApplied,
const G4VParticleChange particleChangeProduced 
)
privatevirtual

Definition at line 192 of file GB01BOptrChangeCrossSection.cc.

199{
200 G4BOptnChangeCrossSection* operation = fChangeCrossSectionOperations[callingProcess];
201 if ( operation == occurenceOperationApplied ) operation->SetInteractionOccured();
202}

Member Data Documentation

◆ fChangeCrossSectionOperations

std::map< const G4BiasingProcessInterface*, G4BOptnChangeCrossSection* > GB01BOptrChangeCrossSection::fChangeCrossSectionOperations
private

Definition at line 99 of file GB01BOptrChangeCrossSection.hh.

◆ fSetup

G4bool GB01BOptrChangeCrossSection::fSetup
private

Definition at line 100 of file GB01BOptrChangeCrossSection.hh.

◆ fParticleToBias

const G4ParticleDefinition* GB01BOptrChangeCrossSection::fParticleToBias
private

Definition at line 101 of file GB01BOptrChangeCrossSection.hh.


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

Applications | User Support | Publications | Collaboration