#include <G4BOptnForceCommonTruncatedExp.hh>
 | 
|   | G4BOptnForceCommonTruncatedExp (G4String name) | 
|   | 
| virtual  | ~G4BOptnForceCommonTruncatedExp () | 
|   | 
virtual const  
G4VBiasingInteractionLaw *  | ProvideOccurenceBiasingInteractionLaw (const G4BiasingProcessInterface *, G4ForceCondition &) | 
|   | 
| virtual G4double  | ProposeAlongStepLimit (const G4BiasingProcessInterface *) | 
|   | 
| virtual G4GPILSelection  | ProposeGPILSelection (const G4GPILSelection processSelection) | 
|   | 
| virtual G4VParticleChange *  | ApplyFinalStateBiasing (const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &) | 
|   | 
| virtual G4double  | DistanceToApplyOperation (const G4Track *, G4double, G4ForceCondition *) | 
|   | 
| virtual G4VParticleChange *  | GenerateBiasingFinalState (const G4Track *, const G4Step *) | 
|   | 
| G4ILawCommonTruncatedExp *  | GetCommonTruncatedExpLaw () | 
|   | 
| G4ILawForceFreeFlight *  | GetForceFreeFlightLaw () | 
|   | 
| void  | Initialize (const G4Track *) | 
|   | 
| void  | UpdateForStep (const G4Step *) | 
|   | 
| void  | Sample () | 
|   | 
| const G4ThreeVector &  | GetInitialMomentum () const  | 
|   | 
| G4double  | GetMaximumDistance () const  | 
|   | 
| void  | ChooseProcessToApply () | 
|   | 
| const G4VProcess *  | GetProcessToApply () const  | 
|   | 
| void  | AddCrossSection (const G4VProcess *, G4double) | 
|   | 
| size_t  | GetNumberOfSharing () const  | 
|   | 
| void  | SetInteractionOccured (G4bool b) | 
|   | 
| G4bool  | GetInteractionOccured () const  | 
|   | 
|   | G4VBiasingOperation (G4String name) | 
|   | 
| virtual  | ~G4VBiasingOperation () | 
|   | 
| virtual void  | AlongMoveBy (const G4BiasingProcessInterface *, const G4Step *, G4double) | 
|   | 
| const G4String &  | GetName () const  | 
|   | 
| std::size_t  | GetUniqueID () const  | 
|   | 
      
        
          | G4BOptnForceCommonTruncatedExp::G4BOptnForceCommonTruncatedExp  | 
          ( | 
          G4String  | 
          name | ) | 
           | 
        
      
 
Definition at line 34 of file G4BOptnForceCommonTruncatedExp.cc.
   37     fProcessToApply(
nullptr),
 
   38     fInteractionOccured(
false),
 
   39     fMaximumDistance(-1.0)
 
   44   fTotalCrossSection = 0.0;
 
G4VBiasingOperation(G4String name)
 
 
 
 
  
  
      
        
          | G4BOptnForceCommonTruncatedExp::~G4BOptnForceCommonTruncatedExp  | 
          ( | 
           | ) | 
           | 
         
       
   | 
  
virtual   | 
  
 
Definition at line 47 of file G4BOptnForceCommonTruncatedExp.cc.
   49   if ( fCommonTruncatedExpLaw ) 
delete fCommonTruncatedExpLaw;
 
   50   if ( fForceFreeFlightLaw )    
delete fForceFreeFlightLaw;
 
 
 
 
Definition at line 116 of file G4BOptnForceCommonTruncatedExp.cc.
  118   fTotalCrossSection      += crossSection;
 
  119   fCrossSections[process]  = crossSection;
 
  120   fNumberOfSharing         = fCrossSections.size();
 
 
 
 
Implements G4VBiasingOperation.
Definition at line 76 of file G4BOptnForceCommonTruncatedExp.cc.
   83       forceFinalState = 
true;
 
   85       return &fDummyParticleChange; 
 
   87   if ( fInteractionOccured )
 
   89       forceFinalState = 
true;
 
   91       return &fDummyParticleChange;
 
   97   if ( processGPIL <= step->GetStepLength() )
 
  103       forceFinalState     = 
false;
 
  104       fInteractionOccured = 
true;
 
  109       forceFinalState = 
true;
 
  111       return &fDummyParticleChange; 
 
virtual void Initialize(const G4Track &track)
 
G4double GetPostStepGPIL() const 
 
G4VProcess * GetWrappedProcess() const 
 
G4double GetAlongStepGPIL() const 
 
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
 
 
 
 
      
        
          | void G4BOptnForceCommonTruncatedExp::ChooseProcessToApply  | 
          ( | 
           | ) | 
           | 
        
      
 
Definition at line 167 of file G4BOptnForceCommonTruncatedExp.cc.
  171   for ( std::map< const G4VProcess*, G4double>::const_iterator it = fCrossSections.begin();
 
  172     it != fCrossSections.end();
 
  175       sigmaSelect += (*it).second;
 
  176       if ( sigmaRand <= sigmaSelect )
 
  178       fProcessToApply = (*it).first;
 
 
 
 
  
  
      
        
          | const G4ThreeVector& G4BOptnForceCommonTruncatedExp::GetInitialMomentum  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
 
  
  
      
        
          | G4bool G4BOptnForceCommonTruncatedExp::GetInteractionOccured  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
 
  
  
      
        
          | G4double G4BOptnForceCommonTruncatedExp::GetMaximumDistance  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
 
  
  
      
        
          | size_t G4BOptnForceCommonTruncatedExp::GetNumberOfSharing  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
 
  
  
      
        
          | const G4VProcess* G4BOptnForceCommonTruncatedExp::GetProcessToApply  | 
          ( | 
           | ) | 
           const | 
         
       
   | 
  
inline   | 
  
 
 
      
        
          | void G4BOptnForceCommonTruncatedExp::Initialize  | 
          ( | 
          const G4Track *  | 
          track | ) | 
           | 
        
      
 
Definition at line 124 of file G4BOptnForceCommonTruncatedExp.cc.
  126   fCrossSections.clear();
 
  127   fTotalCrossSection  = 0.0;
 
  128   fNumberOfSharing    = 0;
 
  130   fInteractionOccured = 
false;
 
  135                   GetNavigatorForTracking()->
 
  136                   GetGlobalToLocalTransform()).TransformPoint(track->
GetPosition());
 
  138                   GetNavigatorForTracking()->
 
  140   fMaximumDistance = currentSolid->
DistanceToOut(localPosition, localDirection);
 
  141   if ( fMaximumDistance <= 
DBL_MIN )  fMaximumDistance = 0.0;
 
const G4ThreeVector & GetPosition() const 
 
G4VSolid * GetSolid() const 
 
static G4TransportationManager * GetTransportationManager()
 
G4ThreeVector GetMomentum() const 
 
const G4ThreeVector & GetMomentumDirection() const 
 
G4LogicalVolume * GetLogicalVolume() const 
 
G4VPhysicalVolume * GetVolume() const 
 
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
 
void SetMaximumDistance(G4double d)
 
 
 
 
      
        
          | void G4BOptnForceCommonTruncatedExp::Sample  | 
          ( | 
           | ) | 
           | 
        
      
 
Definition at line 158 of file G4BOptnForceCommonTruncatedExp.cc.
  161   fCommonTruncatedExpLaw->
Sample();
 
void SetSelectedProcessXSfraction(G4double fXS)
 
void ChooseProcessToApply()
 
void SetForceCrossSection(G4double xs)
 
 
 
 
  
  
      
        
          | void G4BOptnForceCommonTruncatedExp::SetInteractionOccured  | 
          ( | 
          G4bool  | 
          b | ) | 
           | 
         
       
   | 
  
inline   | 
  
 
 
      
        
          | void G4BOptnForceCommonTruncatedExp::UpdateForStep  | 
          ( | 
          const G4Step *  | 
          step | ) | 
           | 
        
      
 
Definition at line 146 of file G4BOptnForceCommonTruncatedExp.cc.
  148   fCrossSections.clear();
 
  149   fTotalCrossSection  = 0.0;
 
  150   fNumberOfSharing    = 0;
 
G4double GetStepLength() const 
 
G4double UpdateForStep(G4double truePathLength)
 
G4double GetMaximumDistance() const 
 
 
 
 
The documentation for this class was generated from the following files: