Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4BiasingProcessInterface Class Reference

#include <G4BiasingProcessInterface.hh>

Inheritance diagram for G4BiasingProcessInterface:
G4VProcess

Public Member Functions

 G4BiasingProcessInterface (G4String name="biasWrapper(0)")
 
 G4BiasingProcessInterface (G4VProcess *wrappedProcess, G4bool wrappedIsAtRest, G4bool wrappedIsAlongStep, G4bool wrappedIsPostStep, G4String useThisName="")
 
 ~G4BiasingProcessInterface ()
 
G4VProcessGetWrappedProcess () const
 
G4VBiasingOperatorGetCurrentBiasingOperator () const
 
G4VBiasingOperatorGetPreviousBiasingOperator () const
 
G4VBiasingOperationGetCurrentOccurenceBiasingOperation () const
 
G4VBiasingOperationGetPreviousOccurenceBiasingOperation () const
 
G4VBiasingOperationGetCurrentFinalStateBiasingOperation () const
 
G4VBiasingOperationGetPreviousFinalStateBiasingOperation () const
 
G4VBiasingOperationGetCurrentNonPhysicsBiasingOperation () const
 
G4VBiasingOperationGetPreviousNonPhysicsBiasingOperation () const
 
G4bool GetIsFirstPostStepGPILInterface (G4bool physOnly=true) const
 
G4bool GetIsLastPostStepGPILInterface (G4bool physOnly=true) const
 
G4bool GetIsFirstPostStepDoItInterface (G4bool physOnly=true) const
 
G4bool GetIsLastPostStepDoItInterface (G4bool physOnly=true) const
 
G4bool IsFirstPostStepGPILInterface (G4bool physOnly=true) const
 
G4bool IsLastPostStepGPILInterface (G4bool physOnly=true) const
 
G4bool IsFirstPostStepDoItInterface (G4bool physOnly=true) const
 
G4bool IsLastPostStepDoItInterface (G4bool physOnly=true) const
 
G4bool GetWrappedProcessIsAtRest () const
 
G4bool GetWrappedProcessIsAlong () const
 
G4bool GetWrappedProcessIsPost () const
 
G4double GetPreviousStepSize () const
 
G4double GetCurrentMinimumStep () const
 
G4double GetProposedSafety () const
 
void SetProposedSafety (G4double sft)
 
G4double GetPostStepGPIL () const
 
G4double GetAlongStepGPIL () const
 
void StartTracking (G4Track *track)
 
void EndTracking ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &step)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &step)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &pd)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &pd)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &pd)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *pd, const G4String &s, G4bool f)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *pd, const G4String &s, G4bool f)
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &pd)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &pd)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
const G4VProcessGetMasterProcess () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 68 of file G4BiasingProcessInterface.hh.

Constructor & Destructor Documentation

G4BiasingProcessInterface::G4BiasingProcessInterface ( G4String  name = "biasWrapper(0)")

Definition at line 45 of file G4BiasingProcessInterface.cc.

References G4Cache< VALTYPE >::Put().

46  : G4VProcess( name ),
47  fCurrentBiasingOperator ( 0 ),
48  fPreviousBiasingOperator( 0 ),
49  fWrappedProcess ( 0 ),
50  fIsPhysicsBasedBiasing ( false ),
51  fWrappedProcessIsAtRest( false ),
52  fWrappedProcessIsAlong ( false ),
53  fWrappedProcessIsPost ( false ),
54  fWrappedProcessInteractionLength( -1.0 ),
55  fBiasingInteractionLaw ( 0 ),
56  fPhysicalInteractionLaw( 0 ),
57  fOccurenceBiasingParticleChange( 0 ),
58  fIamFirstGPIL ( false )
59 {
60  for (G4int i = 0 ; i < 8 ; i++) fFirstLastFlags[i] = false;
61  fResetInteractionLaws.Put( true );
62  fCommonStart.Put(true);
63  fCommonEnd.Put(true);
64 }
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
int G4int
Definition: G4Types.hh:78
void Put(const value_type &val) const
Definition: G4Cache.hh:257
G4BiasingProcessInterface::G4BiasingProcessInterface ( G4VProcess wrappedProcess,
G4bool  wrappedIsAtRest,
G4bool  wrappedIsAlongStep,
G4bool  wrappedIsPostStep,
G4String  useThisName = "" 
)

Definition at line 67 of file G4BiasingProcessInterface.cc.

References G4VProcess::GetProcessName(), G4VProcess::GetProcessSubType(), and G4VProcess::SetProcessSubType().

70  : G4VProcess( useThisName != "" ? useThisName : "biasWrapper("+wrappedProcess->GetProcessName()+")",
71  wrappedProcess->GetProcessType()),
72  fCurrentBiasingOperator ( 0 ),
73  fPreviousBiasingOperator( 0 ),
74  fWrappedProcess ( wrappedProcess ),
75  fIsPhysicsBasedBiasing ( true ),
76  fWrappedProcessIsAtRest( wrappedIsAtRest ),
77  fWrappedProcessIsAlong ( wrappedIsAlongStep ),
78  fWrappedProcessIsPost ( wrappedIsPostStep ),
79  fWrappedProcessInteractionLength( -1.0 ),
80  fBiasingInteractionLaw ( 0 ),
81  fPhysicalInteractionLaw( 0 ),
82  fOccurenceBiasingParticleChange( 0 ),
83  fIamFirstGPIL ( false )
84 {
85  for (G4int i = 0 ; i < 8 ; i++) fFirstLastFlags[i] = false;
86 
87  SetProcessSubType(fWrappedProcess->GetProcessSubType());
88 
89  // -- create physical interaction law:
90  fPhysicalInteractionLaw = new G4InteractionLawPhysical("PhysicalInteractionLawFor("+GetProcessName()+")");
91  // -- instantiate particle change wrapper for occurence biaising:
92  fOccurenceBiasingParticleChange = new G4ParticleChangeForOccurenceBiasing("biasingPCfor"+GetProcessName());
93  fParticleChange = new G4ParticleChange();
94  // -- instantiate a "do nothing" particle change:
95  fDummyParticleChange = new G4ParticleChangeForNothing();
96 }
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
int G4int
Definition: G4Types.hh:78
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
G4BiasingProcessInterface::~G4BiasingProcessInterface ( )

Definition at line 100 of file G4BiasingProcessInterface.cc.

101 {
102  if ( fPhysicalInteractionLaw != 0 ) delete fPhysicalInteractionLaw;
103  if ( fOccurenceBiasingParticleChange ) delete fOccurenceBiasingParticleChange;
104  if ( fDummyParticleChange ) delete fDummyParticleChange;
105 }

Member Function Documentation

G4VParticleChange * G4BiasingProcessInterface::AlongStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 513 of file G4BiasingProcessInterface.cc.

References G4VBiasingOperation::AlongMoveBy(), G4VProcess::AlongStepDoIt(), G4InteractionLawPhysical::ComputeNonInteractionProbabilityAt(), G4endl, G4Exception(), G4VBiasingInteractionLaw::GetName(), G4Step::GetStepLength(), G4Track::GetTrackStatus(), G4ParticleChangeForNothing::Initialize(), JustWarning, G4VParticleChange::ProposeTrackStatus(), G4ParticleChangeForOccurenceBiasing::SetOccurenceWeightForNonInteraction(), and G4ParticleChangeForOccurenceBiasing::SetWrappedParticleChange().

515 {
516  // ---------------------------------------
517  // -- case outside of volume with biasing:
518  // ---------------------------------------
519  if ( fCurrentBiasingOperator == 0 )
520  {
521  if ( fWrappedProcessIsAlong ) return fWrappedProcess->AlongStepDoIt(track, step);
522  else
523  {
524  fDummyParticleChange->Initialize( track );
525  return fDummyParticleChange;
526  }
527  }
528 
529  // -----------------------------------
530  // -- case inside volume with biasing:
531  // -----------------------------------
532  if ( fWrappedProcessIsAlong ) fOccurenceBiasingParticleChange->SetWrappedParticleChange( fWrappedProcess->AlongStepDoIt(track, step) );
533  else
534  {
535  fOccurenceBiasingParticleChange->SetWrappedParticleChange ( 0 );
536  fOccurenceBiasingParticleChange->ProposeTrackStatus( track.GetTrackStatus() );
537  }
538  G4double weightForNonInteraction (1.0);
539  if ( fBiasingInteractionLaw != 0 )
540  {
541  weightForNonInteraction =
542  fPhysicalInteractionLaw->ComputeNonInteractionProbabilityAt(step.GetStepLength()) /
543  fBiasingInteractionLaw ->ComputeNonInteractionProbabilityAt(step.GetStepLength());
544 
545  fOccurenceBiasingOperation->AlongMoveBy( this, &step, weightForNonInteraction );
546 
547  if ( weightForNonInteraction <= 0. )
548  {
550  ed << " Negative non interaction weight : w_NI = " << weightForNonInteraction <<
551  " p_NI(phys) = " << fPhysicalInteractionLaw->ComputeNonInteractionProbabilityAt(step.GetStepLength()) <<
552  " p_NI(bias) = " << fBiasingInteractionLaw ->ComputeNonInteractionProbabilityAt(step.GetStepLength()) <<
553  " step length = " << step.GetStepLength() <<
554  " biasing interaction law = `" << fBiasingInteractionLaw->GetName() << "'" << G4endl;
555  G4Exception(" G4BiasingProcessInterface::AlongStepDoIt(...)",
556  "BIAS.GEN.04",
557  JustWarning,
558  ed);
559  }
560 
561  }
562 
563  fOccurenceBiasingParticleChange->SetOccurenceWeightForNonInteraction( weightForNonInteraction );
564 
565  return fOccurenceBiasingParticleChange;
566 
567 }
virtual void Initialize(const G4Track &track)
virtual G4double ComputeNonInteractionProbabilityAt(G4double length) const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetStepLength() const
G4TrackStatus GetTrackStatus() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
virtual void AlongMoveBy(const G4BiasingProcessInterface *, const G4Step *, G4double)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4double G4BiasingProcessInterface::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 436 of file G4BiasingProcessInterface.cc.

References G4VProcess::AlongStepGetPhysicalInteractionLength(), DBL_MAX, NotCandidateForSelection, G4VBiasingOperation::ProposeAlongStepLimit(), and G4VBiasingOperation::ProposeGPILSelection().

441 {
442  // -- for helper methods:
443  fCurrentMinimumStep = currentMinimumStep;
444  fProposedSafety = proposedSafety;
445 
446 
447  // -- initialization default case:
448  fWrappedProcessAlongStepGPIL = DBL_MAX;
449  *selection = NotCandidateForSelection;
450  // ---------------------------------------
451  // -- case outside of volume with biasing:
452  // ---------------------------------------
453  if ( fCurrentBiasingOperator == 0 )
454  {
455  if ( fWrappedProcessIsAlong ) fWrappedProcessAlongStepGPIL =
456  fWrappedProcess->AlongStepGetPhysicalInteractionLength(track,
457  previousStepSize,
458  currentMinimumStep,
459  proposedSafety,
460  selection);
461  return fWrappedProcessAlongStepGPIL;
462  }
463 
464  // --------------------------------------------------------------------
465  // -- non-physics based biasing: no along operation expected (for now):
466  // --------------------------------------------------------------------
467  if ( !fIsPhysicsBasedBiasing ) return fWrappedProcessAlongStepGPIL;
468 
469  // ----------------------
470  // -- physics-based case:
471  // ----------------------
472  if ( fOccurenceBiasingOperation == 0 )
473  {
474  if ( fWrappedProcessIsAlong ) fWrappedProcessAlongStepGPIL =
475  fWrappedProcess->AlongStepGetPhysicalInteractionLength(track,
476  previousStepSize,
477  currentMinimumStep,
478  proposedSafety,
479  selection);
480  return fWrappedProcessAlongStepGPIL;
481  }
482 
483 
484  // ----------------------------------------------------------
485  // -- From here we have an valid occurence biasing operation:
486  // ----------------------------------------------------------
487  // -- Give operation opportunity to shorten step proposed by physics process:
488  fBiasingAlongStepGPIL = fOccurenceBiasingOperation->ProposeAlongStepLimit( this );
489  G4double minimumStep = fBiasingAlongStepGPIL < currentMinimumStep ? fBiasingAlongStepGPIL : currentMinimumStep ;
490  // -- wrapped process is called with minimum step ( <= currentMinimumStep passed ) : an along process can not
491  // -- have its operation stretched over what it expects:
492  if ( fWrappedProcessIsAlong )
493  {
494  fWrappedProcessAlongStepGPIL = fWrappedProcess->AlongStepGetPhysicalInteractionLength(track,
495  previousStepSize,
496  minimumStep,
497  proposedSafety,
498  selection);
499  fWrappedProcessGPILSelection = *selection;
500  fBiasingGPILSelection = fOccurenceBiasingOperation->ProposeGPILSelection( fWrappedProcessGPILSelection );
501  }
502  else
503  {
504  fBiasingGPILSelection = fOccurenceBiasingOperation->ProposeGPILSelection( NotCandidateForSelection );
505  fWrappedProcessAlongStepGPIL = fBiasingAlongStepGPIL;
506  }
507 
508  *selection = fBiasingGPILSelection;
509  return fWrappedProcessAlongStepGPIL;
510 
511 }
virtual G4double ProposeAlongStepLimit(const G4BiasingProcessInterface *)
virtual G4GPILSelection ProposeGPILSelection(const G4GPILSelection wrappedProcessSelection)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4VParticleChange * G4BiasingProcessInterface::AtRestDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 575 of file G4BiasingProcessInterface.cc.

References G4VProcess::AtRestDoIt().

577 {
578  return fWrappedProcess->AtRestDoIt(track, step);
579 }
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0
G4double G4BiasingProcessInterface::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 570 of file G4BiasingProcessInterface.cc.

References G4VProcess::AtRestGetPhysicalInteractionLength().

572 {
573  return fWrappedProcess->AtRestGetPhysicalInteractionLength(track, condition);
574 }
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)=0
void G4BiasingProcessInterface::BuildPhysicsTable ( const G4ParticleDefinition pd)
virtual

Reimplemented from G4VProcess.

Definition at line 603 of file G4BiasingProcessInterface.cc.

References G4VProcess::BuildPhysicsTable(), and G4VBiasingOperator::GetBiasingOperators().

604 {
605  // -- Inform existing operators about start of the run.
606  // -- IMPORTANT : as PreparePhysicsTable(...) has been called first for all processes,
607  // -- the first/last flags and G4BiasingProcessInterface vector of processes have
608  // -- been properly setup.
609  if ( fIamFirstGPIL )
610  {
611  for ( size_t optr = 0 ; optr < ( G4VBiasingOperator::GetBiasingOperators() ).size() ; optr ++)
612  ( G4VBiasingOperator::GetBiasingOperators() )[optr]->StartRun( );
613  }
614 
615  if ( fWrappedProcess != 0 )
616  {
617  fWrappedProcess->BuildPhysicsTable(pd);
618  }
619 }
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:210
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators()
void G4BiasingProcessInterface::BuildWorkerPhysicsTable ( const G4ParticleDefinition pd)
virtual

Reimplemented from G4VProcess.

Definition at line 654 of file G4BiasingProcessInterface.cc.

References G4VProcess::BuildWorkerPhysicsTable(), and G4VBiasingOperator::GetBiasingOperators().

655 {
656  // -- Inform existing operators about start of the run.
657  // -- IMPORTANT : as PreparePhysicsTable(...) has been called first for all processes,
658  // -- the first/last flags and G4BiasingProcessInterface vector of processes have
659  // -- been properly setup.
660  if ( fIamFirstGPIL )
661  {
662  for ( size_t optr = 0 ; optr < ( G4VBiasingOperator::GetBiasingOperators() ).size() ; optr ++)
663  ( G4VBiasingOperator::GetBiasingOperators() )[optr]->StartRun( );
664  }
665 
666  if ( fWrappedProcess != 0 )
667  {
668  fWrappedProcess->BuildWorkerPhysicsTable(pd);
669  }
670 }
virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition &part)
Definition: G4VProcess.cc:202
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators()
void G4BiasingProcessInterface::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 139 of file G4BiasingProcessInterface.cc.

References G4VProcess::EndTracking(), G4VBiasingOperator::ExitingBiasing(), fKillTrackAndSecondaries, fStopAndKill, G4Cache< VALTYPE >::Get(), G4VBiasingOperator::GetBiasingOperators(), G4BiasingTrackDataStore::GetBiasingTrackData(), G4BiasingTrackDataStore::GetInstance(), G4Step::GetSecondary(), G4Track::GetStep(), G4Track::GetTrackStatus(), and G4Cache< VALTYPE >::Put().

140 {
141  if ( fIsPhysicsBasedBiasing ) fWrappedProcess->EndTracking();
142  if ( fCurrentBiasingOperator) fCurrentBiasingOperator->ExitingBiasing( fCurrentTrack, this );
143  fCurrentBiasingOperator = 0;
144  fPreviousBiasingOperator = 0;
145  fBiasingInteractionLaw = 0;
146 
147 
148  // -- !! this part might have to be improved : could be time consuming
149  // -- !! and assumes all tracks are killed during tracking, which is
150  // -- !! not true : stacking operations may kill tracks
151  if ( fCommonEnd.Get() )
152  {
153  fCommonEnd.Put( false );// = false;
154  fCommonStart.Put( true );// = true;
155 
156  for ( size_t optr = 0 ; optr < ( G4VBiasingOperator::GetBiasingOperators() ).size() ; optr ++)
158 
159  if ( ( fCurrentTrack->GetTrackStatus() == fStopAndKill ) || ( fCurrentTrack->GetTrackStatus() == fKillTrackAndSecondaries ) )
160  {
162  if ( biasingData ) delete biasingData; // -- this also deregisters the biasing data from the track data store
163  if ( fCurrentTrack->GetTrackStatus() == fKillTrackAndSecondaries )
164  {
165  const G4TrackVector* secondaries = fCurrentTrack->GetStep()->GetSecondary();
166  for ( size_t i2nd = 0 ; i2nd < secondaries->size() ; i2nd++ )
167  {
168  biasingData = G4BiasingTrackDataStore::GetInstance()->GetBiasingTrackData( (*secondaries)[i2nd] );
169  if ( biasingData ) delete biasingData;
170  }
171  }
172  }
173  }
174 }
static G4BiasingTrackDataStore * GetInstance()
value_type & Get() const
Definition: G4Cache.hh:253
G4TrackStatus GetTrackStatus() const
const G4Step * GetStep() const
G4BiasingTrackData * GetBiasingTrackData(const G4Track *track)
std::vector< G4Track * > G4TrackVector
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators()
virtual void EndTracking()
Definition: G4VProcess.cc:113
void ExitingBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
void Put(const value_type &val) const
Definition: G4Cache.hh:257
const G4TrackVector * GetSecondary() const
G4double G4BiasingProcessInterface::GetAlongStepGPIL ( ) const
inline

Definition at line 136 of file G4BiasingProcessInterface.hh.

Referenced by G4BOptnForceCommonTruncatedExp::DenyProcessPostStepDoIt().

136 { return fWrappedProcessAlongStepGPIL; }
G4VBiasingOperator* G4BiasingProcessInterface::GetCurrentBiasingOperator ( ) const
inline

Definition at line 91 of file G4BiasingProcessInterface.hh.

91 {return fCurrentBiasingOperator; }
G4VBiasingOperation* G4BiasingProcessInterface::GetCurrentFinalStateBiasingOperation ( ) const
inline

Definition at line 96 of file G4BiasingProcessInterface.hh.

96 { return fFinalStateBiasingOperation; }
G4double G4BiasingProcessInterface::GetCurrentMinimumStep ( ) const
inline

Definition at line 131 of file G4BiasingProcessInterface.hh.

131 { return fCurrentMinimumStep;}
G4VBiasingOperation* G4BiasingProcessInterface::GetCurrentNonPhysicsBiasingOperation ( ) const
inline

Definition at line 98 of file G4BiasingProcessInterface.hh.

98 { return fNonPhysicsBiasingOperation; }
G4VBiasingOperation* G4BiasingProcessInterface::GetCurrentOccurenceBiasingOperation ( ) const
inline

Definition at line 94 of file G4BiasingProcessInterface.hh.

94 { return fOccurenceBiasingOperation; }
G4bool G4BiasingProcessInterface::GetIsFirstPostStepDoItInterface ( G4bool  physOnly = true) const

Definition at line 696 of file G4BiasingProcessInterface.cc.

697 {
698  G4int iPhys = ( physOnly ) ? 1 : 0;
699  return fFirstLastFlags[IdxFirstLast( 1, 0, iPhys)];
700 }
int G4int
Definition: G4Types.hh:78
G4bool G4BiasingProcessInterface::GetIsFirstPostStepGPILInterface ( G4bool  physOnly = true) const

Definition at line 686 of file G4BiasingProcessInterface.cc.

687 {
688  G4int iPhys = ( physOnly ) ? 1 : 0;
689  return fFirstLastFlags[IdxFirstLast( 1, 1, iPhys)];
690 }
int G4int
Definition: G4Types.hh:78
G4bool G4BiasingProcessInterface::GetIsLastPostStepDoItInterface ( G4bool  physOnly = true) const

Definition at line 701 of file G4BiasingProcessInterface.cc.

702 {
703  G4int iPhys = ( physOnly ) ? 1 : 0;
704  return fFirstLastFlags[IdxFirstLast( 0, 0, iPhys)];
705 }
int G4int
Definition: G4Types.hh:78
G4bool G4BiasingProcessInterface::GetIsLastPostStepGPILInterface ( G4bool  physOnly = true) const

Definition at line 691 of file G4BiasingProcessInterface.cc.

692 {
693  G4int iPhys = ( physOnly ) ? 1 : 0;
694  return fFirstLastFlags[IdxFirstLast( 0, 1, iPhys)];
695 }
int G4int
Definition: G4Types.hh:78
G4double G4BiasingProcessInterface::GetPostStepGPIL ( ) const
inline

Definition at line 135 of file G4BiasingProcessInterface.hh.

Referenced by G4BOptnForceCommonTruncatedExp::DenyProcessPostStepDoIt().

135 { return fBiasingPostStepGPIL; }
G4VBiasingOperator* G4BiasingProcessInterface::GetPreviousBiasingOperator ( ) const
inline

Definition at line 92 of file G4BiasingProcessInterface.hh.

92 {return fPreviousBiasingOperator; }
G4VBiasingOperation* G4BiasingProcessInterface::GetPreviousFinalStateBiasingOperation ( ) const
inline

Definition at line 97 of file G4BiasingProcessInterface.hh.

97 { return fPreviousFinalStateBiasingOperation; }
G4VBiasingOperation* G4BiasingProcessInterface::GetPreviousNonPhysicsBiasingOperation ( ) const
inline

Definition at line 99 of file G4BiasingProcessInterface.hh.

99 { return fPreviousNonPhysicsBiasingOperation; }
G4VBiasingOperation* G4BiasingProcessInterface::GetPreviousOccurenceBiasingOperation ( ) const
inline

Definition at line 95 of file G4BiasingProcessInterface.hh.

95 { return fPreviousOccurenceBiasingOperation; }
G4double G4BiasingProcessInterface::GetPreviousStepSize ( ) const
inline

Definition at line 130 of file G4BiasingProcessInterface.hh.

130 { return fPreviousStepSize;}
const G4ProcessManager * G4BiasingProcessInterface::GetProcessManager ( )
virtual

Reimplemented from G4VProcess.

Definition at line 649 of file G4BiasingProcessInterface.cc.

References G4VProcess::GetProcessManager().

650 {
651  if ( fWrappedProcess != 0 ) return fWrappedProcess->GetProcessManager();
652  else return G4VProcess::GetProcessManager();
653 }
virtual const G4ProcessManager * GetProcessManager()
Definition: G4VProcess.hh:514
G4double G4BiasingProcessInterface::GetProposedSafety ( ) const
inline

Definition at line 132 of file G4BiasingProcessInterface.hh.

132 { return fProposedSafety;}
G4VProcess* G4BiasingProcessInterface::GetWrappedProcess ( ) const
inline
G4bool G4BiasingProcessInterface::GetWrappedProcessIsAlong ( ) const
inline

Definition at line 125 of file G4BiasingProcessInterface.hh.

125 { return fWrappedProcessIsAlong; }
G4bool G4BiasingProcessInterface::GetWrappedProcessIsAtRest ( ) const
inline

Definition at line 124 of file G4BiasingProcessInterface.hh.

124 { return fWrappedProcessIsAtRest; }
G4bool G4BiasingProcessInterface::GetWrappedProcessIsPost ( ) const
inline

Definition at line 126 of file G4BiasingProcessInterface.hh.

126 { return fWrappedProcessIsPost; }
G4bool G4BiasingProcessInterface::IsApplicable ( const G4ParticleDefinition pd)
virtual

Reimplemented from G4VProcess.

Definition at line 582 of file G4BiasingProcessInterface.cc.

References G4VProcess::IsApplicable().

583 {
584  if ( fWrappedProcess != 0 ) return fWrappedProcess->IsApplicable(pd);
585  else return true;
586 }
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:205
G4bool G4BiasingProcessInterface::IsFirstPostStepDoItInterface ( G4bool  physOnly = true) const

Definition at line 752 of file G4BiasingProcessInterface.cc.

References G4ProcessManager::GetPostStepProcessVector(), GetWrappedProcess(), G4ProcessVector::size(), and typeDoIt.

753 {
754  G4bool isFirst = true;
755  const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeDoIt);
756  G4int thisIdx(-1);
757  for (G4int i = 0; i < pv->size(); i++ ) if ( (*pv)(i) == this ) { thisIdx = i; break; }
758  for ( size_t i = 0; i < fCoInterfaces->size(); i++ )
759  {
760  if ( ( (*fCoInterfaces)[i]->GetWrappedProcess() != 0 ) || !physOnly )
761  {
762  G4int thatIdx(-1);
763  for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (*fCoInterfaces)[i] ) { thatIdx = j; break; }
764  if ( thisIdx > thatIdx )
765  {
766  isFirst = false;
767  break;
768  }
769  }
770  }
771  return isFirst;
772 }
int G4int
Definition: G4Types.hh:78
G4VProcess * GetWrappedProcess() const
bool G4bool
Definition: G4Types.hh:79
G4int size() const
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4bool G4BiasingProcessInterface::IsFirstPostStepGPILInterface ( G4bool  physOnly = true) const

Definition at line 708 of file G4BiasingProcessInterface.cc.

References G4ProcessManager::GetPostStepProcessVector(), GetWrappedProcess(), G4ProcessVector::size(), and typeGPIL.

709 {
710  G4bool isFirst = true;
711  const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeGPIL);
712  G4int thisIdx(-1);
713  for (G4int i = 0; i < pv->size(); i++ ) if ( (*pv)(i) == this ) { thisIdx = i; break; }
714  for ( size_t i = 0; i < fCoInterfaces->size(); i++ )
715  {
716  if ( ( (*fCoInterfaces)[i]->GetWrappedProcess() != 0 ) || !physOnly )
717  {
718  G4int thatIdx(-1);
719  for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (*fCoInterfaces)[i] ) { thatIdx = j; break; }
720  if ( thisIdx > thatIdx )
721  {
722  isFirst = false;
723  break;
724  }
725  }
726  }
727  return isFirst;
728 }
int G4int
Definition: G4Types.hh:78
G4VProcess * GetWrappedProcess() const
bool G4bool
Definition: G4Types.hh:79
G4int size() const
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4bool G4BiasingProcessInterface::IsLastPostStepDoItInterface ( G4bool  physOnly = true) const

Definition at line 774 of file G4BiasingProcessInterface.cc.

References G4ProcessManager::GetPostStepProcessVector(), GetWrappedProcess(), G4ProcessVector::size(), and typeDoIt.

775 {
776  G4bool isLast = true;
777  const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeDoIt);
778  G4int thisIdx(-1);
779  for (G4int i = 0; i < pv->size(); i++ ) if ( (*pv)(i) == this ) { thisIdx = i; break; }
780  for ( size_t i = 0; i < fCoInterfaces->size(); i++ )
781  {
782  if ( ( (*fCoInterfaces)[i]->GetWrappedProcess() != 0 ) || !physOnly )
783  {
784  G4int thatIdx(-1);
785  for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (*fCoInterfaces)[i] ) { thatIdx = j; break; }
786  if ( thisIdx < thatIdx )
787  {
788  isLast = false;
789  break;
790  }
791  }
792  }
793  return isLast;
794 }
int G4int
Definition: G4Types.hh:78
G4VProcess * GetWrappedProcess() const
bool G4bool
Definition: G4Types.hh:79
G4int size() const
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4bool G4BiasingProcessInterface::IsLastPostStepGPILInterface ( G4bool  physOnly = true) const

Definition at line 730 of file G4BiasingProcessInterface.cc.

References G4ProcessManager::GetPostStepProcessVector(), GetWrappedProcess(), G4ProcessVector::size(), and typeGPIL.

731 {
732  G4bool isLast = true;
733  const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeGPIL);
734  G4int thisIdx(-1);
735  for (G4int i = 0; i < pv->size(); i++ ) if ( (*pv)(i) == this ) { thisIdx = i; break; }
736  for ( size_t i = 0; i < fCoInterfaces->size(); i++ )
737  {
738  if ( ( (*fCoInterfaces)[i]->GetWrappedProcess() != 0 ) || !physOnly )
739  {
740  G4int thatIdx(-1);
741  for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (*fCoInterfaces)[i] ) { thatIdx = j; break; }
742  if ( thisIdx < thatIdx )
743  {
744  isLast = false;
745  break;
746  }
747  }
748  }
749  return isLast;
750 }
int G4int
Definition: G4Types.hh:78
G4VProcess * GetWrappedProcess() const
bool G4bool
Definition: G4Types.hh:79
G4int size() const
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4VParticleChange * G4BiasingProcessInterface::PostStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 316 of file G4BiasingProcessInterface.cc.

References G4VBiasingOperation::ApplyFinalStateBiasing(), BAC_DenyInteraction, BAC_FinalState, BAC_None, BAC_NonPhysics, G4InteractionLawPhysical::ComputeEffectiveCrossSectionAt(), G4VBiasingInteractionLaw::ComputeEffectiveCrossSectionAt(), G4VBiasingOperation::DenyProcessPostStepDoIt(), G4endl, G4Exception(), G4VBiasingOperation::GenerateBiasingFinalState(), G4VBiasingOperator::GetProposedFinalStateBiasingOperation(), G4Step::GetStepLength(), G4VParticleChange::GetTrackStatus(), G4Track::GetWeight(), G4ParticleChange::Initialize(), G4VBiasingInteractionLaw::IsEffectiveCrossSectionInfinite(), G4VBiasingInteractionLaw::IsSingular(), JustWarning, G4VProcess::PostStepDoIt(), G4VParticleChange::ProposeParentWeight(), G4VParticleChange::ProposeTrackStatus(), G4VBiasingOperator::ReportOperationApplied(), G4ParticleChangeForOccurenceBiasing::SetOccurenceWeightForInteraction(), G4VParticleChange::SetSecondaryWeightByProcess(), G4ParticleChangeForOccurenceBiasing::SetWrappedParticleChange(), and G4ParticleChangeForOccurenceBiasing::StealSecondaries().

318 {
319  // ---------------------------------------
320  // -- case outside of volume with biasing:
321  // ---------------------------------------
322  if ( fCurrentBiasingOperator == 0 ) return fWrappedProcess->PostStepDoIt(track, step);
323 
324  // ----------------------------
325  // -- non-physics biasing case:
326  // ----------------------------
327  if ( !fIsPhysicsBasedBiasing )
328  {
329  G4VParticleChange* particleChange = fNonPhysicsBiasingOperation->GenerateBiasingFinalState( &track, &step );
330  fCurrentBiasingOperator->ReportOperationApplied( this, BAC_NonPhysics, fNonPhysicsBiasingOperation, particleChange );
331  return particleChange;
332  }
333 
334  // -- physics biasing case:
335  // ------------------------
336  // -- It proceeds with the following logic:
337  // -- 1) If an occurence biasing operation exists, it makes the
338  // -- decision about the interaction to happen or not.
339  // -- If the occurence operation refuses the interaction, it
340  // -- can only propose a new track weight.
341  // -- 2) The interaction is produced by the final state biasing
342  // -- operation, if it exists, or by the wrapped process in
343  // -- the other case.
344  // -- Hence 2) happens if 1) decides so
345  // 2) happens if there is no occurence biasing operation
346  if ( fOccurenceBiasingOperation != 0 )
347  {
348  G4double proposedTrackWeight = track.GetWeight();
349  if ( fOccurenceBiasingOperation->DenyProcessPostStepDoIt( this, &track, &step, proposedTrackWeight ) )
350  {
351  fParticleChange->Initialize( track ); // **??** <= might use a light version for particle change here
352  fParticleChange->ProposeParentWeight( proposedTrackWeight );
353  fCurrentBiasingOperator->ReportOperationApplied( this, BAC_DenyInteraction, fOccurenceBiasingOperation, fParticleChange );
354  return fParticleChange;
355  }
356  }
357 
358 
359  // -- usual case with generated final state:
360  G4VParticleChange* finalStateParticleChange;
362  fFinalStateBiasingOperation = fCurrentBiasingOperator->GetProposedFinalStateBiasingOperation( &track, this );
363  if ( fFinalStateBiasingOperation != 0 )
364  {
365  finalStateParticleChange = fFinalStateBiasingOperation->ApplyFinalStateBiasing( this, &track, &step );
366  BAC = BAC_FinalState;
367  }
368  else
369  {
370  finalStateParticleChange = fWrappedProcess->PostStepDoIt(track, step);
371  BAC = BAC_None ;
372  }
373  // -- if no occurence biasing operation, we're done:
374  if ( fOccurenceBiasingOperation == 0 )
375  {
376  fCurrentBiasingOperator->ReportOperationApplied( this, BAC, fFinalStateBiasingOperation, finalStateParticleChange );
377  return finalStateParticleChange;
378  }
379 
380  // -- If occurence biasing, applies on top of final state occurence biasing weight correction:
381  G4double weightForInteraction = 1.0;
382  if ( !fBiasingInteractionLaw->IsSingular() ) weightForInteraction =
383  fPhysicalInteractionLaw->ComputeEffectiveCrossSectionAt(step.GetStepLength()) /
384  fBiasingInteractionLaw ->ComputeEffectiveCrossSectionAt(step.GetStepLength());
385  else
386  {
387  // -- at this point effective XS can only be infinite, if not, there is a logic problem
388  if ( !fBiasingInteractionLaw->IsEffectiveCrossSectionInfinite() )
389  {
391  ed << "Internal inconsistency in cross-section handling. Please report !" << G4endl;
392  G4Exception(" G4BiasingProcessInterface::PostStepDoIt(...)",
393  "BIAS.GEN.02",
394  JustWarning,
395  ed);
396  // -- if XS is infinite, weight is zero (and will stay zero), but we'll do differently.
397  // -- Should foresee in addition something to remember that in case of singular
398  // -- distribution, weight can only be partly calculated
399  }
400  }
401 
402  if ( weightForInteraction <= 0. )
403  {
405  ed << " Negative interaction weight : w_I = "
406  << weightForInteraction <<
407  " XS_I(phys) = " << fBiasingInteractionLaw ->ComputeEffectiveCrossSectionAt(step.GetStepLength()) <<
408  " XS_I(bias) = " << fPhysicalInteractionLaw->ComputeEffectiveCrossSectionAt(step.GetStepLength()) <<
409  " step length = " << step.GetStepLength() <<
410  " Interaction law = `" << fBiasingInteractionLaw << "'" <<
411  G4endl;
412  G4Exception(" G4BiasingProcessInterface::PostStepDoIt(...)",
413  "BIAS.GEN.03",
414  JustWarning,
415  ed);
416 
417  }
418 
419  fCurrentBiasingOperator->ReportOperationApplied( this, BAC,
420  fOccurenceBiasingOperation, weightForInteraction,
421  fFinalStateBiasingOperation, finalStateParticleChange );
422 
423  fOccurenceBiasingParticleChange->SetOccurenceWeightForInteraction( weightForInteraction );
424  fOccurenceBiasingParticleChange->SetSecondaryWeightByProcess( true );
425  fOccurenceBiasingParticleChange->SetWrappedParticleChange( finalStateParticleChange );
426  fOccurenceBiasingParticleChange->ProposeTrackStatus( finalStateParticleChange->GetTrackStatus() );
427  fOccurenceBiasingParticleChange->StealSecondaries(); // -- this also makes weightForInteraction applied to secondaries stolen
428 
429  // -- finish:
430  return fOccurenceBiasingParticleChange;
431 
432 }
G4VBiasingOperation * GetProposedFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetStepLength() const
virtual G4bool DenyProcessPostStepDoIt(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4double &)
void ProposeParentWeight(G4double finalWeight)
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *)=0
virtual G4double ComputeEffectiveCrossSectionAt(G4double length) const =0
G4BiasingAppliedCase
void SetSecondaryWeightByProcess(G4bool)
virtual G4double ComputeEffectiveCrossSectionAt(G4double length) const
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *)=0
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4bool IsSingular() const
virtual void Initialize(const G4Track &)
G4double GetWeight() const
#define G4endl
Definition: G4ios.hh:61
virtual G4bool IsEffectiveCrossSectionInfinite() const
G4TrackStatus GetTrackStatus() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void ReportOperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
G4double G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 178 of file G4BiasingProcessInterface.cc.

References condition(), DBL_MAX, G4VBiasingOperation::DistanceToApplyOperation(), G4VBiasingOperator::ExitingBiasing(), fGeomBoundary, G4VBiasingOperator::GetBiasingOperator(), G4VProcess::GetCurrentInteractionLength(), G4Track::GetCurrentStepNumber(), G4VPhysicalVolume::GetLogicalVolume(), G4Step::GetPreStepPoint(), G4VBiasingOperator::GetProposedNonPhysicsBiasingOperation(), G4VBiasingOperator::GetProposedOccurenceBiasingOperation(), G4VBiasingInteractionLaw::GetSampledInteractionLength(), G4Track::GetStep(), G4StepPoint::GetStepStatus(), G4Track::GetVolume(), NotCandidateForSelection, NotForced, G4VProcess::PostStepGetPhysicalInteractionLength(), G4VBiasingOperation::ProposeForceCondition(), G4VBiasingOperation::ProvideOccurenceBiasingInteractionLaw(), G4VProcess::ResetNumberOfInteractionLengthLeft(), and G4InteractionLawPhysical::SetPhysicalCrossSection().

181 {
182  // -- Remember previous operator and proposed operations, if any, and reset:
183  // -------------------------------------------------------------------------
184  // -- remember:
185  fPreviousBiasingOperator = fCurrentBiasingOperator;
186  fPreviousOccurenceBiasingOperation = fOccurenceBiasingOperation;
187  fPreviousFinalStateBiasingOperation = fFinalStateBiasingOperation;
188  fPreviousNonPhysicsBiasingOperation = fNonPhysicsBiasingOperation;
189  fPreviousBiasingInteractionLaw = fBiasingInteractionLaw;
190  // -- reset:
191  fOccurenceBiasingOperation = 0;
192  fFinalStateBiasingOperation = 0;
193  fNonPhysicsBiasingOperation = 0;
194  fBiasingInteractionLaw = 0;
195  // -- Physics PostStep and AlongStep GPIL
196  fWrappedProcessPostStepGPIL = DBL_MAX;
197  fBiasingPostStepGPIL = DBL_MAX;
198  fWrappedProcessInteractionLength = DBL_MAX; // -- inverse of analog cross-section, no biasing counter-part in general
199  fWrappedProcessForceCondition = NotForced;
200  fBiasingForceCondition = NotForced;
201  fWrappedProcessAlongStepGPIL = DBL_MAX;
202  fBiasingAlongStepGPIL = DBL_MAX;
203  fWrappedProcessGPILSelection = NotCandidateForSelection;
204  fBiasingGPILSelection = NotCandidateForSelection;
205  // -- for helper:
206  fPreviousStepSize = previousStepSize;
207 
208 
209  // -- If new volume, get possible new biasing operator:
210  // ----------------------------------------------------
211  // -- [Note : bug with this first step ! Does not work if previous step was concurrently limited with geometry]
212  G4bool firstStepInVolume = ( (track.GetStep()->GetPreStepPoint()->GetStepStatus() == fGeomBoundary) || (track.GetCurrentStepNumber() == 1) );
213  if ( firstStepInVolume ) fCurrentBiasingOperator = G4VBiasingOperator::GetBiasingOperator( track.GetVolume()->GetLogicalVolume() );
214 
215  // -----------------------------------
216  // -- Previous step was under biasing:
217  // -----------------------------------
218  if ( fPreviousBiasingOperator != 0 )
219  {
220  // -- if current step not under same operator, let this operator knows:
221  if ( fPreviousBiasingOperator != fCurrentBiasingOperator ) fPreviousBiasingOperator->ExitingBiasing( &track, this );
222  // -- if biasing does not continue, set process behavior to standard tracking:
223  // -- ... and see [*] below...
224  if ( fCurrentBiasingOperator == 0 )
225  {
226  fResetWrappedProcessInteractionLength = true;
227  ResetForUnbiasedTracking();
228  }
229  }
230 
231  // --------------------------------------------------------------
232  // -- no operator : analog tracking if physics-based, or nothing:
233  // --------------------------------------------------------------
234  if ( fCurrentBiasingOperator == 0 )
235  {
236  if ( fIsPhysicsBasedBiasing )
237  {
238  // -- [*] this wrapped process has just been ResetForUnbiasedTracking(), it has a fresh #int length, and is
239  // -- in a state where it believes it is starting tracking : let it believe hence the previous step
240  // -- length was 0.0, as for a normal first step:
241  if ( fResetWrappedProcessInteractionLength )
242  {
243  fResetWrappedProcessInteractionLength = false;
244  fWrappedProcess->ResetNumberOfInteractionLengthLeft();
245  return fWrappedProcess->PostStepGetPhysicalInteractionLength(track, 0.0, condition);
246  }
247  return fWrappedProcess->PostStepGetPhysicalInteractionLength(track, previousStepSize, condition);
248  }
249  else
250  {
251  *condition = NotForced;
252  return DBL_MAX;
253  }
254  }
255 
256 
257 
258  // --------------------------------------------------
259  // -- An biasing operator exists. Proceed with
260  // -- treating non-physics and physics biasing cases:
261  //---------------------------------------------------
262 
263  // -- non-physics-based biasing case:
264  // ----------------------------------
265  if ( !fIsPhysicsBasedBiasing )
266  {
267  fNonPhysicsBiasingOperation = fCurrentBiasingOperator->GetProposedNonPhysicsBiasingOperation( &track, this );
268  if ( fNonPhysicsBiasingOperation == 0 )
269  {
270  *condition = NotForced;
271  return DBL_MAX;
272  }
273  return fNonPhysicsBiasingOperation->DistanceToApplyOperation(&track, previousStepSize, condition);
274  }
275 
276 
277  // -- Physics based biasing case:
278  // ------------------------------
279  // -- call to underneath physics process PostStepGPIL to update it with current point:
280  fWrappedProcessPostStepGPIL = fWrappedProcess->PostStepGetPhysicalInteractionLength(track, previousStepSize, condition);
281  fWrappedProcessForceCondition = *condition;
282  // -- **! At this point, might have to be careful of previousStepSize being larger than the previous process
283  // -- **! PostStepGPIL proposed: as we disregard this PostStepGPIL value, such situation does happen. Does not
284  // -- **! look to have generated problems for now, but might be fragile.
285 
286  // -- Ask for possible GPIL biasing operation:
287  fOccurenceBiasingOperation = fCurrentBiasingOperator->GetProposedOccurenceBiasingOperation( &track, this );
288 
289  // -- no operation for occurence biasing, analog GPIL returns the wrapped process GPIL and condition values
290  // -- (note that condition was set above):
291  if ( fOccurenceBiasingOperation == 0 ) return fWrappedProcessPostStepGPIL;
292 
293  // -- A valid GPIL biasing operation has been proposed:
294  // -- 0) remember wrapped process will need to be reset on biasing exit, if particle survives:
295  fResetWrappedProcessInteractionLength = true;
296  // -- 1) collect/update process interaction length for reference analog interaction law:
297  fWrappedProcessInteractionLength = fWrappedProcess->GetCurrentInteractionLength();
298  fPhysicalInteractionLaw->SetPhysicalCrossSection( 1.0 / fWrappedProcessInteractionLength );
299  // -- 2) Collect biasing interaction law:
300  // -- The interaction law pointer is collected as a const pointer to the interaction law object.
301  // -- This interaction law will be kept under control of the biasing operation, which is the only
302  // -- entity that will change the state of the biasing interaction law.
303  fBiasingInteractionLaw = fOccurenceBiasingOperation->ProvideOccurenceBiasingInteractionLaw( this );
304  // -- 3) Ask operation to sample the biasing interaction law:
305  fBiasingPostStepGPIL = fBiasingInteractionLaw->GetSampledInteractionLength();
306  fBiasingForceCondition = fOccurenceBiasingOperation->ProposeForceCondition( fWrappedProcessForceCondition );
307 
308  // -- finish
309  *condition = fBiasingForceCondition;
310  return fBiasingPostStepGPIL;
311 
312 }
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
G4StepStatus GetStepStatus() const
G4VBiasingOperation * GetProposedOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
const G4Step * GetStep() const
G4StepPoint * GetPreStepPoint() const
G4int GetCurrentStepNumber() const
bool G4bool
Definition: G4Types.hh:79
virtual G4ForceCondition ProposeForceCondition(const G4ForceCondition wrappedProcessCondition)
G4double GetCurrentInteractionLength() const
Definition: G4VProcess.hh:462
void SetPhysicalCrossSection(G4double crossSection)
static G4VBiasingOperator * GetBiasingOperator(const G4LogicalVolume *)
G4LogicalVolume * GetLogicalVolume() const
G4VPhysicalVolume * GetVolume() const
G4double GetSampledInteractionLength() const
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *)=0
G4VBiasingOperation * GetProposedNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
void ExitingBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
#define DBL_MAX
Definition: templates.hh:83
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *)=0
void G4BiasingProcessInterface::PreparePhysicsTable ( const G4ParticleDefinition pd)
virtual

Reimplemented from G4VProcess.

Definition at line 620 of file G4BiasingProcessInterface.cc.

References G4VProcess::PreparePhysicsTable().

621 {
622  // -- Let process finding its first/last position in the process manager:
623  SetUpFirstLastFlags();
624  if ( fWrappedProcess != 0 )
625  {
626  fWrappedProcess->PreparePhysicsTable(pd);
627  }
628 }
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:217
void G4BiasingProcessInterface::PrepareWorkerPhysicsTable ( const G4ParticleDefinition pd)
virtual

Reimplemented from G4VProcess.

Definition at line 671 of file G4BiasingProcessInterface.cc.

References G4VProcess::PrepareWorkerPhysicsTable().

672 {
673  // -- Let process finding its first/last position in the process manager:
674  SetUpFirstLastFlags();
675  if ( fWrappedProcess != 0 )
676  {
677  fWrappedProcess->PrepareWorkerPhysicsTable(pd);
678  }
679 }
virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.cc:207
void G4BiasingProcessInterface::ResetNumberOfInteractionLengthLeft ( )
virtual

Reimplemented from G4VProcess.

Definition at line 680 of file G4BiasingProcessInterface.cc.

References G4VProcess::ResetNumberOfInteractionLengthLeft().

681 {
682  if ( fWrappedProcess != 0 ) fWrappedProcess->ResetNumberOfInteractionLengthLeft();
683 }
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
G4bool G4BiasingProcessInterface::RetrievePhysicsTable ( const G4ParticleDefinition pd,
const G4String s,
G4bool  f 
)
virtual

Reimplemented from G4VProcess.

Definition at line 634 of file G4BiasingProcessInterface.cc.

References G4VProcess::RetrievePhysicsTable().

635 {
636  if ( fWrappedProcess != 0 ) return fWrappedProcess->RetrievePhysicsTable(pd, s, f);
637  else return false;
638 }
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:236
void G4BiasingProcessInterface::SetMasterProcess ( G4VProcess masterP)
virtual

Reimplemented from G4VProcess.

Definition at line 589 of file G4BiasingProcessInterface.cc.

References G4VProcess::GetMasterProcess(), GetWrappedProcess(), and G4VProcess::SetMasterProcess().

590 {
591  // -- Master for this process:
593  // -- Master for wrapped process:
594  if ( fWrappedProcess != 0 )
595  {
596  const G4BiasingProcessInterface* thisWrapperMaster = (const G4BiasingProcessInterface *)GetMasterProcess();
597  // -- paranoia check:
598  G4VProcess* wrappedMaster = 0;
599  wrappedMaster = thisWrapperMaster->GetWrappedProcess();
600  fWrappedProcess->SetMasterProcess( wrappedMaster );
601  }
602 }
virtual void SetMasterProcess(G4VProcess *masterP)
Definition: G4VProcess.cc:212
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
G4VProcess * GetWrappedProcess() const
void G4BiasingProcessInterface::SetProcessManager ( const G4ProcessManager mgr)
virtual

Reimplemented from G4VProcess.

Definition at line 640 of file G4BiasingProcessInterface.cc.

References G4VProcess::SetProcessManager().

641 {
642  if ( fWrappedProcess != 0 ) fWrappedProcess->SetProcessManager(mgr);
644  (fManagerInterfaceMap[mgr]).push_back(this);
645  fCoInterfaces = &(fManagerInterfaceMap[mgr]);
646  fProcessManager = mgr;
647 }
virtual void SetProcessManager(const G4ProcessManager *)
Definition: G4VProcess.hh:508
void G4BiasingProcessInterface::SetProposedSafety ( G4double  sft)
inline

Definition at line 133 of file G4BiasingProcessInterface.hh.

133 { fProposedSafety = sft;}
void G4BiasingProcessInterface::StartTracking ( G4Track track)
virtual

Reimplemented from G4VProcess.

Definition at line 108 of file G4BiasingProcessInterface.cc.

References G4Cache< VALTYPE >::Get(), G4VBiasingOperator::GetBiasingOperators(), G4Cache< VALTYPE >::Put(), and G4VProcess::StartTracking().

109 {
110  fCurrentTrack = track;
111  if ( fIsPhysicsBasedBiasing ) fWrappedProcess->StartTracking(fCurrentTrack);
112  fCurrentBiasingOperator = 0;
113  fPreviousBiasingOperator = 0;
114  fOccurenceBiasingOperation = 0;
115  fPreviousOccurenceBiasingOperation = 0;
116  fFinalStateBiasingOperation = 0;
117  fPreviousFinalStateBiasingOperation = 0;
118  fNonPhysicsBiasingOperation = 0;
119  fPreviousNonPhysicsBiasingOperation = 0;
120  fBiasingInteractionLaw = 0;
121  fPreviousBiasingInteractionLaw = 0;
122 
123  fPreviousStepSize = -1.0;
124 
125 
126  fResetWrappedProcessInteractionLength = false;
127 
128  if ( fCommonStart.Get() )
129  {
130  fCommonStart.Put( false );// = false;
131  fCommonEnd.Put(true);// = true;
132 
133  for ( size_t optr = 0 ; optr < ( G4VBiasingOperator::GetBiasingOperators() ).size() ; optr ++)
134  ( G4VBiasingOperator::GetBiasingOperators() )[optr]->StartTracking( fCurrentTrack );
135  }
136 }
value_type & Get() const
Definition: G4Cache.hh:253
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:101
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators()
void Put(const value_type &val) const
Definition: G4Cache.hh:257
G4bool G4BiasingProcessInterface::StorePhysicsTable ( const G4ParticleDefinition pd,
const G4String s,
G4bool  f 
)
virtual

Reimplemented from G4VProcess.

Definition at line 629 of file G4BiasingProcessInterface.cc.

References G4VProcess::StorePhysicsTable().

630 {
631  if ( fWrappedProcess != 0 ) return fWrappedProcess->StorePhysicsTable(pd, s, f);
632  else return false;
633 }
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:231

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