Geant4  10.02.p03
GB01BOptrChangeCrossSection Class Reference

#include <GB01BOptrChangeCrossSection.hh>

Inheritance diagram for GB01BOptrChangeCrossSection:
Collaboration diagram for GB01BOptrChangeCrossSection:

Public Member Functions

 GB01BOptrChangeCrossSection (G4String particleToBias, G4String name="ChangeXS")
 
virtual ~GB01BOptrChangeCrossSection ()
 
virtual void StartRun ()
 
- Public Member Functions inherited from G4VBiasingOperator
 G4VBiasingOperator (G4String name)
 
virtual ~G4VBiasingOperator ()
 
virtual void Configure ()
 
virtual void ConfigureForWorker ()
 
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 ()
 

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
 

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 void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void ExitBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 

Detailed Description

Definition at line 53 of file GB01BOptrChangeCrossSection.hh.

Constructor & Destructor Documentation

◆ GB01BOptrChangeCrossSection()

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

Definition at line 38 of file GB01BOptrChangeCrossSection.cc.

40  : G4VBiasingOperator(name),
41  fSetup(true)
42 {
44 
45  if ( fParticleToBias == 0 )
46  {
48  ed << "Particle `" << particleName << "' not found !" << G4endl;
49  G4Exception("GB01BOptrChangeCrossSection(...)",
50  "exGB01.01",
52  ed);
53  }
54 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4ParticleDefinition * fParticleToBias
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4ParticleTable * GetParticleTable()
G4VBiasingOperator(G4String name)
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ ~GB01BOptrChangeCrossSection()

GB01BOptrChangeCrossSection::~GB01BOptrChangeCrossSection ( )
virtual

Definition at line 56 of file GB01BOptrChangeCrossSection.cc.

57 {
58  for ( std::map< const G4BiasingProcessInterface*, G4BOptnChangeCrossSection* >::iterator
59  it = fChangeCrossSectionOperations.begin() ;
60  it != fChangeCrossSectionOperations.end() ;
61  it++ ) delete (*it).second;
62 }
std::map< const G4BiasingProcessInterface *, G4BOptnChangeCrossSection *> fChangeCrossSectionOperations

Member Function Documentation

◆ OperationApplied()

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

Reimplemented from G4VBiasingOperator.

Definition at line 183 of file GB01BOptrChangeCrossSection.cc.

189 {
190  G4BOptnChangeCrossSection* operation = fChangeCrossSectionOperations[callingProcess];
191  if ( operation == occurenceOperationApplied ) operation->SetInteractionOccured();
192 }
std::map< const G4BiasingProcessInterface *, G4BOptnChangeCrossSection *> fChangeCrossSectionOperations
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ProposeFinalStateBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 74 of file GB01BOptrChangeCrossSection.hh.

75  {return 0;}

◆ ProposeNonPhysicsBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 77 of file GB01BOptrChangeCrossSection.hh.

78  {return 0;}
Here is the call graph for this function:

◆ ProposeOccurenceBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 96 of file GB01BOptrChangeCrossSection.cc.

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

◆ StartRun()

void GB01BOptrChangeCrossSection::StartRun ( )
virtual

Reimplemented from G4VBiasingOperator.

Definition at line 65 of file GB01BOptrChangeCrossSection.cc.

66 {
67  // --------------
68  // -- Setup stage:
69  // ---------------
70  // -- Start by collecting processes under biasing, create needed biasing
71  // -- operations and associate these operations to the processes:
72  if ( fSetup )
73  {
74  const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
75  const G4BiasingProcessSharedData* sharedData =
77  if ( sharedData ) // -- sharedData tested, as is can happen a user attaches an operator to a
78  // -- volume but without defined BiasingProcessInterface processes.
79  {
80  for ( size_t i = 0 ; i < (sharedData->GetPhysicsBiasingProcessInterfaces()).size(); i++ )
81  {
82  const G4BiasingProcessInterface* wrapperProcess =
83  (sharedData->GetPhysicsBiasingProcessInterfaces())[i];
84  G4String operationName = "XSchange-" +
85  wrapperProcess->GetWrappedProcess()->GetProcessName();
86  fChangeCrossSectionOperations[wrapperProcess] =
87  new G4BOptnChangeCrossSection(operationName);
88  }
89  }
90  fSetup = false;
91  }
92 }
const std::vector< const G4BiasingProcessInterface *> & GetPhysicsBiasingProcessInterfaces() const
G4ProcessManager * GetProcessManager() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const G4BiasingProcessSharedData * GetSharedData() const
const G4ParticleDefinition * fParticleToBias
G4VProcess * GetWrappedProcess() const
std::map< const G4BiasingProcessInterface *, G4BOptnChangeCrossSection *> fChangeCrossSectionOperations
Here is the call graph for this function:

Member Data Documentation

◆ fChangeCrossSectionOperations

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

Definition at line 98 of file GB01BOptrChangeCrossSection.hh.

◆ fParticleToBias

const G4ParticleDefinition* GB01BOptrChangeCrossSection::fParticleToBias
private

Definition at line 100 of file GB01BOptrChangeCrossSection.hh.

◆ fSetup

G4bool GB01BOptrChangeCrossSection::fSetup
private

Definition at line 99 of file GB01BOptrChangeCrossSection.hh.


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