G4VProcess Class Reference

#include <G4VProcess.hh>

Inheritance diagram for G4VProcess:

G4AdjointProcessEquivalentToDirectProcess G4Cerenkov G4CoupledTransportation G4FastSimulationManagerProcess G4ImportanceProcess G4IVContinuousDiscreteProcess G4IVRestDiscreteProcess G4ParallelWorldProcess G4ParallelWorldScoringProcess G4ScoreSplittingProcess G4StepLimiter G4Transportation G4UserSpecialCuts G4VContinuousDiscreteProcess G4VContinuousProcess G4VDiscreteProcess G4VITProcess G4VRestContinuousDiscreteProcess G4VRestContinuousProcess G4VRestDiscreteProcess G4VRestProcess G4WeightCutOffProcess G4WeightWindowProcess G4WrapperProcess SpecialCuts

Public Member Functions

 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
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
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)
virtual G4bool IsApplicable (const G4ParticleDefinition &)
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
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)
virtual void StartTracking (G4Track *)
virtual void EndTracking ()
virtual void SetProcessManager (const G4ProcessManager *)
virtual const G4ProcessManagerGetProcessManager ()
virtual void ResetNumberOfInteractionLengthLeft ()
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

Static Public Member Functions

static const G4StringGetProcessTypeName (G4ProcessType)

Protected Member Functions

void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
void ClearNumberOfInteractionLengthLeft ()

Protected Attributes

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 75 of file G4VProcess.hh.


Constructor & Destructor Documentation

G4VProcess::G4VProcess ( const G4String aName = "NoName",
G4ProcessType  aType = fNotDefined 
)

Definition at line 51 of file G4VProcess.cc.

References aParticleChange, and pParticleChange.

00052                   : aProcessManager(0),
00053                     pParticleChange(0),
00054                     theNumberOfInteractionLengthLeft(-1.0),
00055                     currentInteractionLength(-1.0),
00056                     theInitialNumberOfInteractionLength(-1.0),
00057                     theProcessName(aName),
00058                     theProcessType(aType),
00059                     theProcessSubType(-1),
00060                     thePILfactor(1.0),
00061                     enableAtRestDoIt(true),
00062                     enableAlongStepDoIt(true),
00063                     enablePostStepDoIt(true),
00064                     verboseLevel(0)
00065 
00066 {
00067   pParticleChange = &aParticleChange;
00068 }

G4VProcess::G4VProcess ( const G4VProcess right  ) 

Definition at line 74 of file G4VProcess.cc.

G4VProcess::~G4VProcess (  )  [virtual]

Definition at line 70 of file G4VProcess.cc.

00071 {
00072 }


Member Function Documentation

virtual G4VParticleChange* G4VProcess::AlongStepDoIt ( const G4Track track,
const G4Step stepData 
) [pure virtual]

Implemented in SpecialCuts, G4ImportanceProcess, G4WeightCutOffProcess, G4WeightWindowProcess, G4AdjointAlongStepWeightCorrection, G4AdjointProcessEquivalentToDirectProcess, G4ContinuousGainOfEnergy, G4ITTransportation, G4VITRestDiscreteProcess, G4VITRestProcess, G4DNABrownianTransportation, G4DNASecondOrderReaction, G4VeLowEnergyLoss, G4ErrorEnergyLoss, G4hImpactIonisation, G4NuclearStopping, G4VEnergyLoss, G4VEnergyLossProcess, G4VMultipleScattering, G4Cerenkov, G4IVContinuousDiscreteProcess, G4IVRestDiscreteProcess, G4VContinuousDiscreteProcess, G4VContinuousProcess, G4VDiscreteProcess, G4VRestContinuousDiscreteProcess, G4VRestContinuousProcess, G4VRestDiscreteProcess, G4VRestProcess, G4WrapperProcess, G4FastSimulationManagerProcess, G4ParallelWorldProcess, G4ParallelWorldScoringProcess, G4ScoreSplittingProcess, G4CoupledTransportation, G4StepLimiter, G4Transportation, and G4UserSpecialCuts.

Referenced by G4WrapperProcess::AlongStepDoIt(), G4AdjointProcessEquivalentToDirectProcess::AlongStepDoIt(), and G4ITStepProcessor::InvokeAlongStepDoItProcs().

virtual G4double G4VProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
) [pure virtual]

Implemented in SpecialCuts, G4ImportanceProcess, G4WeightCutOffProcess, G4WeightWindowProcess, G4AdjointProcessEquivalentToDirectProcess, G4ITTransportation, G4VITRestDiscreteProcess, G4VITRestProcess, G4DNABrownianTransportation, G4DNASecondOrderReaction, G4NuclearStopping, G4VEnergyLossProcess, G4VMultipleScattering, G4Cerenkov, G4IVContinuousDiscreteProcess, G4IVRestDiscreteProcess, G4VContinuousDiscreteProcess, G4VContinuousProcess, G4VDiscreteProcess, G4VRestContinuousDiscreteProcess, G4VRestContinuousProcess, G4VRestDiscreteProcess, G4VRestProcess, G4WrapperProcess, G4FastSimulationManagerProcess, G4ParallelWorldProcess, G4ParallelWorldScoringProcess, G4ScoreSplittingProcess, G4CoupledTransportation, G4StepLimiter, G4Transportation, and G4UserSpecialCuts.

Referenced by AlongStepGPIL().

G4double G4VProcess::AlongStepGPIL ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
) [inline]

Definition at line 450 of file G4VProcess.hh.

References AlongStepGetPhysicalInteractionLength().

00455 {
00456   G4double value
00457    =AlongStepGetPhysicalInteractionLength(track, previousStepSize, currentMinimumStep, proposedSafety, selection);
00458   return value;
00459 }

virtual G4VParticleChange* G4VProcess::AtRestDoIt ( const G4Track track,
const G4Step stepData 
) [pure virtual]

Implemented in SpecialCuts, G4ImportanceProcess, G4WeightCutOffProcess, G4WeightWindowProcess, G4Decay, G4AdjointProcessEquivalentToDirectProcess, G4ITTransportation, G4VITRestDiscreteProcess, G4VITRestProcess, G4DNAMolecularDecay, G4DNASecondOrderReaction, G4eplusPolarizedAnnihilation, G4eplusAnnihilation, G4Cerenkov, G4Scintillation, G4PionMinusNuclearAtRestChips, G4ProtonAntiProtonAtRestChips, G4QCaptureAtRest, G4AntiNeutronAnnihilationAtRest, G4AntiProtonAnnihilationAtRest, G4HadronStoppingProcess, G4KaonMinusAbsorption, G4KaonMinusAbsorptionAtRest, G4MuonMinusCaptureAtRest, G4NeutronCaptureAtRest, G4PiMinusAbsorptionAtRest, G4PionMinusAbsorptionAtRest, G4IVContinuousDiscreteProcess, G4IVRestDiscreteProcess, G4VContinuousDiscreteProcess, G4VContinuousProcess, G4VDiscreteProcess, G4VRestContinuousDiscreteProcess, G4VRestContinuousProcess, G4VRestDiscreteProcess, G4VRestProcess, G4WrapperProcess, G4FastSimulationManagerProcess, G4ParallelWorldProcess, G4ParallelWorldScoringProcess, G4ScoreSplittingProcess, G4CoupledTransportation, G4StepLimiter, G4Transportation, and G4UserSpecialCuts.

Referenced by G4WrapperProcess::AtRestDoIt(), G4AdjointProcessEquivalentToDirectProcess::AtRestDoIt(), and G4ITStepProcessor::InvokeAtRestDoItProcs().

virtual G4double G4VProcess::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
) [pure virtual]

Implemented in SpecialCuts, G4ImportanceProcess, G4WeightCutOffProcess, G4WeightWindowProcess, G4Decay, G4AdjointProcessEquivalentToDirectProcess, G4ITTransportation, G4VITRestDiscreteProcess, G4VITRestProcess, G4DNAMolecularDecay, G4DNASecondOrderReaction, G4eplusPolarizedAnnihilation, G4eplusAnnihilation, G4Cerenkov, G4PionMinusNuclearAtRestChips, G4ProtonAntiProtonAtRestChips, G4AntiNeutronAnnihilationAtRest, G4AntiProtonAnnihilationAtRest, G4HadronStoppingProcess, G4KaonMinusAbsorption, G4NeutronCaptureAtRest, G4PionMinusAbsorptionAtRest, G4IVContinuousDiscreteProcess, G4IVRestDiscreteProcess, G4VContinuousDiscreteProcess, G4VContinuousProcess, G4VDiscreteProcess, G4VRestContinuousDiscreteProcess, G4VRestContinuousProcess, G4VRestDiscreteProcess, G4VRestProcess, G4WrapperProcess, G4FastSimulationManagerProcess, G4ParallelWorldProcess, G4ParallelWorldScoringProcess, G4ScoreSplittingProcess, G4CoupledTransportation, G4StepLimiter, G4Transportation, and G4UserSpecialCuts.

Referenced by G4WrapperProcess::AtRestGetPhysicalInteractionLength(), G4AdjointProcessEquivalentToDirectProcess::AtRestGetPhysicalInteractionLength(), and AtRestGPIL().

G4double G4VProcess::AtRestGPIL ( const G4Track track,
G4ForceCondition condition 
) [inline]

Definition at line 461 of file G4VProcess.hh.

References AtRestGetPhysicalInteractionLength(), and thePILfactor.

Referenced by G4ITStepProcessor::GetAtRestIL().

00463 {
00464   G4double value
00465    =AtRestGetPhysicalInteractionLength(track, condition);
00466   return thePILfactor*value;
00467 }

virtual void G4VProcess::BuildPhysicsTable ( const G4ParticleDefinition  )  [inline, virtual]

Reimplemented in G4Decay, G4UnknownDecay, G4AdjointAlongStepWeightCorrection, G4AdjointProcessEquivalentToDirectProcess, G4ContinuousGainOfEnergy, G4VAdjointReverseReaction, G4ITTransportation, G4VITProcess, G4DNABrownianTransportation, G4DNASecondOrderReaction, G4AnnihiToMuPair, G4GammaConversionToMuons, G4VLowEnergyDiscretePhotonProcess, G4hImpactIonisation, G4eplusPolarizedAnnihilation, G4ePolarizedIonisation, G4PolarizedCompton, G4VEmProcess, G4VEnergyLossProcess, G4VMultipleScattering, G4SynchrotronRadiation, G4VXTRenergyLoss, G4HadronicProcess, G4PionMinusNuclearAtRestChips, G4ProtonAntiProtonAtRestChips, G4QInelastic, G4QNGamma, G4ChargeExchangeProcess, G4RadioactiveDecay, G4AntiNeutronAnnihilationAtRest, G4AntiProtonAnnihilationAtRest, G4HadronStoppingProcess, G4KaonMinusAbsorption, G4KaonMinusAbsorptionAtRest, G4MuonMinusCaptureAtRest, G4NeutronCaptureAtRest, G4PiMinusAbsorptionAtRest, G4PionMinusAbsorptionAtRest, and G4WrapperProcess.

Definition at line 210 of file G4VProcess.hh.

Referenced by G4VUserPhysicsList::BuildIntegralPhysicsTable(), G4WrapperProcess::BuildPhysicsTable(), and G4AdjointProcessEquivalentToDirectProcess::BuildPhysicsTable().

00210 {}

void G4VProcess::ClearNumberOfInteractionLengthLeft (  )  [inline, protected]

Reimplemented in G4VITProcess.

Definition at line 418 of file G4VProcess.hh.

References theInitialNumberOfInteractionLength, and theNumberOfInteractionLengthLeft.

Referenced by G4VContinuousDiscreteProcess::AlongStepDoIt(), G4IVContinuousDiscreteProcess::AlongStepDoIt(), G4VRestProcess::AtRestDoIt(), G4VRestDiscreteProcess::AtRestDoIt(), G4VRestContinuousProcess::AtRestDoIt(), G4VRestContinuousDiscreteProcess::AtRestDoIt(), G4IVRestDiscreteProcess::AtRestDoIt(), G4UnknownDecay::DecayIt(), G4RadioactiveDecay::DecayIt(), G4Decay::DecayIt(), G4Decay::EndTracking(), G4WHadronElasticProcess::PostStepDoIt(), G4VRestDiscreteProcess::PostStepDoIt(), G4VRestContinuousDiscreteProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VDiscreteProcess::PostStepDoIt(), G4VContinuousDiscreteProcess::PostStepDoIt(), G4VAdjointReverseReaction::PostStepDoIt(), G4TransitionRadiation::PostStepDoIt(), G4IVRestDiscreteProcess::PostStepDoIt(), G4IVContinuousDiscreteProcess::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), and G4HadronElasticProcess::PostStepDoIt().

00419 {
00420   theInitialNumberOfInteractionLength = -1.0; 
00421   theNumberOfInteractionLengthLeft =  -1.0;
00422 }

void G4VProcess::DumpInfo (  )  const [virtual]

Definition at line 206 of file G4VProcess.cc.

References G4cout, G4endl, GetProcessTypeName(), theProcessName, theProcessSubType, and theProcessType.

Referenced by G4ProcessTable::DumpInfo(), and G4ProcessManagerMessenger::SetNewValue().

00207 {
00208   G4cout << "Process Name " << theProcessName ;
00209   G4cout << " : Type[" << GetProcessTypeName(theProcessType) << "]";
00210   G4cout << " : SubType[" << theProcessSubType << "]"<< G4endl;
00211 }

void G4VProcess::EndTracking (  )  [virtual]

Reimplemented in G4Decay, G4AdjointProcessEquivalentToDirectProcess, G4WrapperProcess, G4FastSimulationManagerProcess, and G4CoupledTransportation.

Definition at line 137 of file G4VProcess.cc.

References currentInteractionLength, G4cout, G4endl, theInitialNumberOfInteractionLength, theNumberOfInteractionLengthLeft, theProcessName, and verboseLevel.

Referenced by G4WrapperProcess::EndTracking(), and G4AdjointProcessEquivalentToDirectProcess::EndTracking().

00138 {
00139 #ifdef G4VERBOSE
00140   if (verboseLevel>2) {
00141     G4cout << "G4VProcess::EndTracking() [" << theProcessName << "]" <<G4endl;
00142   }
00143 #endif
00144   theNumberOfInteractionLengthLeft = -1.0;
00145   currentInteractionLength = -1.0;
00146   theInitialNumberOfInteractionLength=-1.0;
00147 }

G4double G4VProcess::GetCurrentInteractionLength (  )  const [inline]

Definition at line 433 of file G4VProcess.hh.

References currentInteractionLength.

00434 {
00435   return currentInteractionLength;
00436 }

G4double G4VProcess::GetNumberOfInteractionLengthLeft (  )  const [inline]

Definition at line 424 of file G4VProcess.hh.

References theNumberOfInteractionLengthLeft.

00425 {
00426   return theNumberOfInteractionLengthLeft;
00427 }

const G4String & G4VProcess::GetPhysicsTableFileName ( const G4ParticleDefinition ,
const G4String directory,
const G4String tableName,
G4bool  ascii = false 
)

Definition at line 214 of file G4VProcess.cc.

References G4ParticleDefinition::GetParticleName(), thePhysicsTableFileName, and theProcessName.

Referenced by G4VEmProcess::RetrievePhysicsTable(), G4VMultipleScattering::StorePhysicsTable(), and G4VEmProcess::StorePhysicsTable().

00218 {
00219   G4String thePhysicsTableFileExt;
00220   if (ascii) thePhysicsTableFileExt = ".asc";
00221   else       thePhysicsTableFileExt = ".dat";
00222 
00223   thePhysicsTableFileName = directory + "/";
00224   thePhysicsTableFileName += tableName + "." +  theProcessName + ".";
00225   thePhysicsTableFileName += particle->GetParticleName() + thePhysicsTableFileExt;
00226   
00227   return thePhysicsTableFileName;
00228 }

G4double G4VProcess::GetPILfactor (  )  const [inline]

Definition at line 445 of file G4VProcess.hh.

References thePILfactor.

00446 {
00447   return thePILfactor;
00448 }

const G4ProcessManager * G4VProcess::GetProcessManager (  )  [inline, virtual]

Reimplemented in G4WrapperProcess.

Definition at line 485 of file G4VProcess.hh.

References aProcessManager.

Referenced by G4WrapperProcess::GetProcessManager().

00486 {
00487   return  aProcessManager; 
00488 }

const G4String & G4VProcess::GetProcessName (  )  const [inline]

Definition at line 379 of file G4VProcess.hh.

References theProcessName.

Referenced by G4VEnergyLossProcess::ActivateForcedInteraction(), G4VEmProcess::ActivateForcedInteraction(), G4VEnergyLossProcess::ActivateSecondaryBiasing(), G4VEmProcess::ActivateSecondaryBiasing(), G4VEnergyLossProcess::AddCollaborativeProcess(), G4ProcessManager::AddProcess(), G4SteppingVerbose::AlongStepDoItAllDone(), G4SteppingVerbose::AlongStepDoItOneByOne(), G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4IVContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4SteppingVerbose::AtRestDoItInvoked(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestProcess::AtRestGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4ProtonAntiProtonAtRestChips::AtRestGetPhysicalInteractionLength(), G4PionMinusNuclearAtRestChips::AtRestGetPhysicalInteractionLength(), G4PionMinusAbsorptionAtRest::AtRestGetPhysicalInteractionLength(), G4NeutronCaptureAtRest::AtRestGetPhysicalInteractionLength(), G4KaonMinusAbsorption::AtRestGetPhysicalInteractionLength(), G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4AntiProtonAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), G4AntiNeutronAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), G4HadronicProcess::BiasCrossSectionByFactor(), G4QInelasticCHIPSBuilder::Build(), G4QInelasticCHIPS_HPBuilder::Build(), G4MiscQGSCBuilder::Build(), G4MiscCHIPSBuilder::Build(), G4HyperonCHIPSBuilder::Build(), G4VXTRenergyLoss::BuildAngleForEnergyBank(), G4VEnergyLossProcess::BuildDEDXTable(), G4VUserPhysicsList::BuildIntegralPhysicsTable(), G4VEnergyLossProcess::BuildLambdaTable(), G4VMultipleScattering::BuildPhysicsTable(), G4VEnergyLossProcess::BuildPhysicsTable(), G4VEmProcess::BuildPhysicsTable(), G4LossTableManager::BuildPhysicsTable(), G4DNABrownianTransportation::BuildPhysicsTable(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4ProcessManager::CheckOrderingParameters(), G4HadronicProcess::CheckResult(), G4StackChecker::ClassifyNewTrack(), G4HadronHElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronDElasticPhysics::ConstructProcess(), G4RichTrajectoryPoint::CreateAttValues(), G4RichTrajectory::CreateAttValues(), G4WHadronElasticProcess::Description(), G4HadronElasticProcess::Description(), G4SteppingVerbose::DPSLAlongStep(), G4SteppingVerbose::DPSLPostStep(), G4HadronicProcessStore::Dump(), G4HadronicProcess::DumpState(), G4AntiNeutronAnnihilationAtRest::G4AntiNeutronAnnihilationAtRest(), G4AntiProtonAnnihilationAtRest::G4AntiProtonAnnihilationAtRest(), G4Cerenkov::G4Cerenkov(), G4ErrorEnergyLoss::G4ErrorEnergyLoss(), G4ErrorTrackLengthTarget::G4ErrorTrackLengthTarget(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4ImportanceProcess::G4ImportanceProcess(), G4KaonMinusAbsorption::G4KaonMinusAbsorption(), G4KaonMinusAbsorptionAtRest::G4KaonMinusAbsorptionAtRest(), G4NeutronCaptureAtRest::G4NeutronCaptureAtRest(), G4OpAbsorption::G4OpAbsorption(), G4OpBoundaryProcess::G4OpBoundaryProcess(), G4OpMieHG::G4OpMieHG(), G4OpRayleigh::G4OpRayleigh(), G4OpWLS::G4OpWLS(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess(), G4PiMinusAbsorptionAtRest::G4PiMinusAbsorptionAtRest(), G4PionMinusAbsorptionAtRest::G4PionMinusAbsorptionAtRest(), G4QAtomicElectronScattering::G4QAtomicElectronScattering(), G4QCaptureAtRest::G4QCaptureAtRest(), G4QCoherentChargeExchange::G4QCoherentChargeExchange(), G4QDiffraction::G4QDiffraction(), G4QDiscProcessMixer::G4QDiscProcessMixer(), G4QElastic::G4QElastic(), G4QInelastic::G4QInelastic(), G4QIonIonElastic::G4QIonIonElastic(), G4QLowEnergy::G4QLowEnergy(), G4QNGamma::G4QNGamma(), G4Scintillation::G4Scintillation(), G4ScoreSplittingProcess::G4ScoreSplittingProcess(), G4StepLimiter::G4StepLimiter(), G4UserSpecialCuts::G4UserSpecialCuts(), G4VLowEnergyDiscretePhotonProcess::G4VLowEnergyDiscretePhotonProcess(), G4WeightCutOffProcess::G4WeightCutOffProcess(), G4WeightWindowProcess::G4WeightWindowProcess(), G4WeightWindowProcess::GetName(), G4ProcTblElement::GetProcessName(), G4ProcessManager::GetProcessVectorIndex(), G4hhIonisation::InitialiseEnergyLossProcess(), G4ProcessTable::Insert(), G4ErrorPropagator::MakeOneStep(), MaxTimeCuts::MaxTimeCuts(), MinEkineCuts::MinEkineCuts(), G4VEmProcess::PostStepDoIt(), G4SteppingVerbose::PostStepDoItAllDone(), G4SteppingVerbose::PostStepDoItOneByOne(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4EmConfigurator::PrepareModels(), G4VMultipleScattering::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), G4VEmProcess::PreparePhysicsTable(), G4LossTableManager::PreparePhysicsTable(), G4HadronicProcessStore::PrintHtml(), G4VMultipleScattering::PrintInfoDefinition(), G4VEnergyLossProcess::PrintInfoDefinition(), G4VEmProcess::PrintInfoDefinition(), G4SynchrotronRadiation::PrintInfoDefinition(), G4hImpactIonisation::PrintInfoDefinition(), G4GammaConversionToMuons::PrintInfoDefinition(), G4AnnihiToMuPair::PrintInfoDefinition(), G4LossTableManager::Register(), G4LossTableManager::RegisterExtraParticle(), G4HadronicProcess::RegisterMe(), G4WrapperProcess::RegisterProcess(), G4PhysicsListHelper::RegisterProcess(), G4ProcessTable::Remove(), G4ProcessPlacer::RemoveProcess(), G4ProcessManager::RemoveProcess(), G4VEnergyLossProcess::RetrievePhysicsTable(), G4VEmProcess::RetrievePhysicsTable(), G4VEnergyLossProcess::SetCrossSectionBiasingFactor(), G4VEmProcess::SetCrossSectionBiasingFactor(), G4VEnergyLossProcess::SetInverseRangeTable(), G4VEnergyLossProcess::SetLambdaTable(), G4ProcessTableMessenger::SetNewValue(), G4ProcessTable::SetProcessActivation(), G4ProcessManager::SetProcessOrdering(), G4ProcessManager::SetProcessOrderingToFirst(), G4ProcessManager::SetProcessOrderingToLast(), G4ProcessManager::SetProcessOrderingToSecond(), G4VEnergyLossProcess::SetRangeTableForLoss(), G4VEnergyLossProcess::SetSecondaryRangeTable(), G4VEnergyLossProcess::SetSubLambdaTable(), G4FastSimulationManagerProcess::SetWorldVolume(), G4SteppingVerbose::ShowStep(), SpecialCuts::SpecialCuts(), G4SteppingVerbose::StepInfo(), G4VMultipleScattering::StorePhysicsTable(), G4VEnergyLossProcess::StorePhysicsTable(), G4VEmProcess::StorePhysicsTable(), G4ScoreSplittingProcess::Verbose(), G4ParallelWorldScoringProcess::Verbose(), G4SteppingVerbose::VerboseTrack(), G4VEmProcess::~G4VEmProcess(), G4VEnergyLossProcess::~G4VEnergyLossProcess(), and G4VMultipleScattering::~G4VMultipleScattering().

00380 {
00381   return theProcessName;
00382 }

G4int G4VProcess::GetProcessSubType (  )  const [inline]

Definition at line 397 of file G4VProcess.hh.

References theProcessSubType.

Referenced by G4VMultipleScattering::BuildPhysicsTable(), G4DNABrownianTransportation::BuildPhysicsTable(), G4HadronicProcessStore::FindProcess(), G4HadronicProcessStore::GetCrossSectionPerAtom(), G4HadronicProcessStore::GetCrossSectionPerVolume(), G4VMultipleScattering::PrintInfoDefinition(), G4VEnergyLossProcess::PrintInfoDefinition(), G4VEmProcess::PrintInfoDefinition(), G4GammaConversionToMuons::PrintInfoDefinition(), G4AnnihiToMuPair::PrintInfoDefinition(), G4PhysicsListHelper::RegisterProcess(), and G4ElectronIonPair::ResidualeChargePostStep().

00398 {
00399   return theProcessSubType;
00400 }

G4ProcessType G4VProcess::GetProcessType (  )  const [inline]

Definition at line 385 of file G4VProcess.hh.

References theProcessType.

Referenced by G4RichTrajectoryPoint::CreateAttValues(), G4RichTrajectory::CreateAttValues(), G4AdjointProcessEquivalentToDirectProcess::G4AdjointProcessEquivalentToDirectProcess(), G4WrapperProcess::RegisterProcess(), G4PhysicsListHelper::RegisterProcess(), and G4ProcessTableMessenger::SetNewValue().

00386 {
00387   return theProcessType;
00388 }

const G4String & G4VProcess::GetProcessTypeName ( G4ProcessType   )  [static]

Definition at line 150 of file G4VProcess.cc.

References fDecay, fElectromagnetic, fGeneral, fHadronic, fNotDefined, fOptical, fParameterisation, fPhotolepton_hadron, fTransportation, and fUserDefined.

Referenced by G4RichTrajectoryPoint::CreateAttValues(), G4RichTrajectory::CreateAttValues(), DumpInfo(), G4ProcessManager::DumpInfo(), G4ProcessTableMessenger::G4ProcessTableMessenger(), and G4ProcessTableMessenger::GetCurrentValue().

00151 {
00152   static G4String typeNotDefined= "NotDefined";
00153   static G4String typeTransportation = "Transportation";
00154   static G4String typeElectromagnetic = "Electromagnetic";
00155   static G4String typeOptical = "Optical";
00156   static G4String typeHadronic = "Hadronic";
00157   static G4String typePhotolepton_hadron = "Photolepton_hadron";
00158   static G4String typeDecay = "Decay";
00159   static G4String typeGeneral = "General";
00160   static G4String typeParameterisation = "Parameterisation";
00161   static G4String typeUserDefined = "UserDefined";
00162   static G4String noType = "------";   // Do not modify this !!!!
00163 
00164   if (aType ==   fNotDefined) {
00165     return  typeNotDefined;
00166   } else if  (aType ==   fTransportation ) {
00167     return typeTransportation;
00168   } else if  (aType ==   fElectromagnetic ) {
00169     return typeElectromagnetic;
00170   } else if  (aType ==   fOptical ) {
00171     return typeOptical;
00172   } else if  (aType ==   fHadronic ) {
00173     return typeHadronic;
00174   } else if  (aType ==   fPhotolepton_hadron ) {
00175     return typePhotolepton_hadron;
00176   } else if  (aType ==   fDecay ) {
00177     return typeDecay;
00178   } else if  (aType ==   fGeneral ) {
00179     return typeGeneral;
00180   } else if  (aType ==   fParameterisation ) {
00181     return typeParameterisation;
00182   } else if  (aType ==   fUserDefined ) {
00183     return typeUserDefined;
00184   } else {
00185     return noType;  
00186   }
00187 }

G4double G4VProcess::GetTotalNumberOfInteractionLengthTraversed (  )  const [inline]

Definition at line 429 of file G4VProcess.hh.

References theInitialNumberOfInteractionLength, and theNumberOfInteractionLengthLeft.

G4int G4VProcess::GetVerboseLevel (  )  const [inline]

Reimplemented in G4Decay, G4UnknownDecay, G4ITTransportation, G4RadioactiveDecay, G4CoupledTransportation, and G4Transportation.

Definition at line 413 of file G4VProcess.hh.

References verboseLevel.

Referenced by G4ProcessTable::DumpInfo().

00414 {
00415   return  verboseLevel;
00416 }

G4bool G4VProcess::isAlongStepDoItIsEnabled (  )  const [inline]

Definition at line 497 of file G4VProcess.hh.

References enableAlongStepDoIt.

Referenced by G4ProcessManager::CheckOrderingParameters().

00498 {
00499   return enableAlongStepDoIt;
00500 }

virtual G4bool G4VProcess::IsApplicable ( const G4ParticleDefinition  )  [inline, virtual]

Reimplemented in G4Decay, G4UnknownDecay, G4AdjointhMultipleScattering, G4AdjointProcessEquivalentToDirectProcess, G4DNAAttachment, G4DNAChargeDecrease, G4DNAChargeIncrease, G4DNAElastic, G4DNAElectronSolvatation, G4DNAExcitation, G4DNAIonisation, G4DNAMolecularDecay, G4DNAVibExcitation, G4AnnihiToMuPair, G4eeToHadrons, G4GammaConversionToMuons, G4hBremsstrahlung, G4hhIonisation, G4hPairProduction, G4mplIonisation, G4MuElecElastic, G4MuElecInelastic, G4RayleighScattering, G4VLowEnergyDiscretePhotonProcess, G4ErrorEnergyLoss, G4MuBremsstrahlung, G4MuIonisation, G4MuMultipleScattering, G4MuPairProduction, G4hImpactIonisation, G4eplusPolarizedAnnihilation, G4ePolarizedIonisation, G4PolarizedCompton, G4PolarizedGammaConversion, G4PolarizedPhotoElectricEffect, G4alphaIonisation, G4ComptonScattering, G4CoulombScattering, G4eBremsstrahlung, G4eIonisation, G4eMultipleScattering, G4eplusAnnihilation, G4GammaConversion, G4hIonisation, G4hMultipleScattering, G4ionIonisation, G4NuclearStopping, G4PhotoElectricEffect, G4VEmProcess, G4VEnergyLossProcess, G4VMultipleScattering, G4Cerenkov, G4Scintillation, G4SynchrotronRadiation, G4SynchrotronRadiationInMat, G4TransitionRadiation, G4VTransitionRadiation, G4VXTRenergyLoss, G4HadronInelasticProcess, G4PionMinusNuclearAtRestChips, G4ProtonAntiProtonAtRestChips, G4QAtomicElectronScattering, G4QCaptureAtRest, G4QCoherentChargeExchange, G4QDiffraction, G4QDiscProcessMixer, G4QElastic, G4QInelastic, G4QIonIonElastic, G4QLowEnergy, G4QNGamma, G4QSynchRad, G4ChargeExchangeProcess, G4RadioactiveDecay, G4HadronCaptureProcess, G4HadronFissionProcess, G4MuonNuclearProcess, G4AntiNeutronAnnihilationAtRest, G4AntiProtonAnnihilationAtRest, G4HadronicAbsorptionBertini, G4HadronicAbsorptionFritiof, G4HadronStoppingProcess, G4KaonMinusAbsorption, G4KaonMinusAbsorptionAtRest, G4MuonMinusCapture, G4MuonMinusCaptureAtRest, G4NeutronCaptureAtRest, G4PiMinusAbsorptionAtRest, G4PionMinusAbsorptionAtRest, G4WrapperProcess, G4OpAbsorption, G4OpBoundaryProcess, G4OpMieHG, G4OpRayleigh, G4OpWLS, and G4NeutronKiller.

Definition at line 205 of file G4VProcess.hh.

Referenced by G4ProcessManager::AddProcess(), G4WrapperProcess::IsApplicable(), and G4AdjointProcessEquivalentToDirectProcess::IsApplicable().

00205 {return true;}

G4bool G4VProcess::isAtRestDoItIsEnabled (  )  const [inline]

Definition at line 491 of file G4VProcess.hh.

References enableAtRestDoIt.

Referenced by G4ProcessManager::CheckOrderingParameters().

00492 {
00493   return enableAtRestDoIt;
00494 }

G4bool G4VProcess::isPostStepDoItIsEnabled (  )  const [inline]

Definition at line 503 of file G4VProcess.hh.

References enablePostStepDoIt.

Referenced by G4ProcessManager::CheckOrderingParameters().

00504 {
00505   return enablePostStepDoIt;
00506 }

G4int G4VProcess::operator!= ( const G4VProcess right  )  const

Definition at line 201 of file G4VProcess.cc.

00202 {
00203   return (this !=  &right);
00204 }

G4int G4VProcess::operator== ( const G4VProcess right  )  const

Definition at line 196 of file G4VProcess.cc.

00197 {
00198   return (this == &right);
00199 }

virtual G4VParticleChange* G4VProcess::PostStepDoIt ( const G4Track track,
const G4Step stepData 
) [pure virtual]

Implemented in G4ErrorTrackLengthTarget, G4VErrorLimitProcess, SpecialCuts, G4ImportanceProcess, G4WeightCutOffProcess, G4WeightWindowProcess, G4Decay, G4UnknownDecay, G4AdjointProcessEquivalentToDirectProcess, G4VAdjointReverseReaction, G4ITTransportation, G4VITRestDiscreteProcess, G4VITRestProcess, G4DNABrownianTransportation, G4DNASecondOrderReaction, G4AnnihiToMuPair, G4GammaConversionToMuons, G4VeLowEnergyLoss, G4hImpactIonisation, G4hRDEnergyLoss, G4VEmProcess, G4VEnergyLoss, G4VEnergyLossProcess, G4VMultipleScattering, G4Cerenkov, G4ForwardXrayTR, G4Scintillation, G4SynchrotronRadiation, G4SynchrotronRadiationInMat, G4TransitionRadiation, G4VTransitionRadiation, G4VXTRenergyLoss, G4HadronicProcess, G4QAtomicElectronScattering, G4QCoherentChargeExchange, G4QDiffraction, G4QDiscProcessMixer, G4QElastic, G4QInelastic, G4QIonIonElastic, G4QLowEnergy, G4QNGamma, G4QSynchRad, G4HadronElasticProcess, G4WHadronElasticProcess, G4IVContinuousDiscreteProcess, G4IVRestDiscreteProcess, G4VContinuousDiscreteProcess, G4VContinuousProcess, G4VDiscreteProcess, G4VRestContinuousDiscreteProcess, G4VRestContinuousProcess, G4VRestDiscreteProcess, G4VRestProcess, G4WrapperProcess, G4OpAbsorption, G4OpBoundaryProcess, G4OpMieHG, G4OpRayleigh, G4OpWLS, G4FastSimulationManagerProcess, G4ParallelWorldProcess, G4ParallelWorldScoringProcess, G4ScoreSplittingProcess, G4CoupledTransportation, G4NeutronKiller, G4StepLimiter, G4Transportation, and G4UserSpecialCuts.

Referenced by G4ITStepProcessor::InvokePSDIP(), G4WrapperProcess::PostStepDoIt(), and G4AdjointProcessEquivalentToDirectProcess::PostStepDoIt().

virtual G4double G4VProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
) [pure virtual]

Implemented in G4ErrorMagFieldLimitProcess, G4ErrorStepLengthLimitProcess, G4ErrorTrackLengthTarget, G4VErrorLimitProcess, MaxTimeCuts, MinEkineCuts, SpecialCuts, G4ImportanceProcess, G4WeightCutOffProcess, G4WeightWindowProcess, G4Decay, G4UnknownDecay, G4AdjointProcessEquivalentToDirectProcess, G4ITTransportation, G4VITRestDiscreteProcess, G4VITRestProcess, G4DNASecondOrderReaction, G4eplusPolarizedAnnihilation, G4ePolarizedIonisation, G4PolarizedCompton, G4VEmProcess, G4VEnergyLossProcess, G4VMultipleScattering, G4Cerenkov, G4QDiscProcessMixer, G4HadronStoppingProcess, G4IVContinuousDiscreteProcess, G4IVRestDiscreteProcess, G4VContinuousDiscreteProcess, G4VContinuousProcess, G4VDiscreteProcess, G4VRestContinuousDiscreteProcess, G4VRestContinuousProcess, G4VRestDiscreteProcess, G4VRestProcess, G4WrapperProcess, G4FastSimulationManagerProcess, G4ParallelWorldProcess, G4ParallelWorldScoringProcess, G4ScoreSplittingProcess, G4CoupledTransportation, G4NeutronKiller, G4StepLimiter, G4Transportation, and G4UserSpecialCuts.

Referenced by G4WrapperProcess::PostStepGetPhysicalInteractionLength(), G4AdjointProcessEquivalentToDirectProcess::PostStepGetPhysicalInteractionLength(), and PostStepGPIL().

G4double G4VProcess::PostStepGPIL ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
) [inline]

Definition at line 469 of file G4VProcess.hh.

References PostStepGetPhysicalInteractionLength(), and thePILfactor.

00472 {
00473   G4double value
00474    =PostStepGetPhysicalInteractionLength(track, previousStepSize, condition);
00475   return thePILfactor*value;
00476 }

virtual void G4VProcess::PreparePhysicsTable ( const G4ParticleDefinition  )  [inline, virtual]

Reimplemented in G4AdjointAlongStepWeightCorrection, G4AdjointProcessEquivalentToDirectProcess, G4ContinuousGainOfEnergy, G4VAdjointReverseReaction, G4eplusPolarizedAnnihilation, G4PolarizedCompton, G4VEmProcess, G4VEnergyLossProcess, G4VMultipleScattering, G4HadronicProcess, G4HadronElasticProcess, G4AntiNeutronAnnihilationAtRest, G4AntiProtonAnnihilationAtRest, G4HadronStoppingProcess, G4KaonMinusAbsorption, G4KaonMinusAbsorptionAtRest, G4MuonMinusCaptureAtRest, G4NeutronCaptureAtRest, G4PiMinusAbsorptionAtRest, G4PionMinusAbsorptionAtRest, and G4WrapperProcess.

Definition at line 217 of file G4VProcess.hh.

Referenced by G4WrapperProcess::PreparePhysicsTable(), and G4AdjointProcessEquivalentToDirectProcess::PreparePhysicsTable().

00217 {}

void G4VProcess::ResetNumberOfInteractionLengthLeft (  )  [virtual]

Reimplemented in G4AdjointProcessEquivalentToDirectProcess, G4VITProcess, and G4WrapperProcess.

Definition at line 92 of file G4VProcess.cc.

References G4UniformRand, theInitialNumberOfInteractionLength, and theNumberOfInteractionLengthLeft.

Referenced by G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4ProtonAntiProtonAtRestChips::AtRestGetPhysicalInteractionLength(), G4PionMinusNuclearAtRestChips::AtRestGetPhysicalInteractionLength(), G4PionMinusAbsorptionAtRest::AtRestGetPhysicalInteractionLength(), G4NeutronCaptureAtRest::AtRestGetPhysicalInteractionLength(), G4KaonMinusAbsorption::AtRestGetPhysicalInteractionLength(), G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4AntiProtonAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), G4AntiNeutronAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4WrapperProcess::ResetNumberOfInteractionLengthLeft(), G4AdjointProcessEquivalentToDirectProcess::ResetNumberOfInteractionLengthLeft(), and G4Decay::StartTracking().

virtual G4bool G4VProcess::RetrievePhysicsTable ( const G4ParticleDefinition ,
const G4String ,
G4bool   
) [inline, virtual]

Reimplemented in G4AdjointProcessEquivalentToDirectProcess, G4VEmProcess, G4VEnergyLossProcess, G4VMultipleScattering, and G4WrapperProcess.

Definition at line 236 of file G4VProcess.hh.

Referenced by G4WrapperProcess::RetrievePhysicsTable(), and G4AdjointProcessEquivalentToDirectProcess::RetrievePhysicsTable().

00237                                                                    {return false;}

void G4VProcess::SetPILfactor ( G4double  value  )  [inline]

Definition at line 438 of file G4VProcess.hh.

References thePILfactor.

00439 {
00440   if (value>0.) {
00441     thePILfactor = value;
00442   }
00443 }

void G4VProcess::SetProcessManager ( const G4ProcessManager  )  [inline, virtual]

Reimplemented in G4WrapperProcess.

Definition at line 479 of file G4VProcess.hh.

References aProcessManager.

Referenced by G4WrapperProcess::SetProcessManager().

00480 {
00481    aProcessManager = procMan; 
00482 }

void G4VProcess::SetProcessSubType ( G4int   )  [inline]

Definition at line 403 of file G4VProcess.hh.

References theProcessSubType.

Referenced by G4alphaIonisation::G4alphaIonisation(), G4AnnihiToMuPair::G4AnnihiToMuPair(), G4AntiNeutronAnnihilationAtRest::G4AntiNeutronAnnihilationAtRest(), G4AntiProtonAnnihilationAtRest::G4AntiProtonAnnihilationAtRest(), G4Cerenkov::G4Cerenkov(), G4ComptonScattering::G4ComptonScattering(), G4CoulombScattering::G4CoulombScattering(), G4CoupledTransportation::G4CoupledTransportation(), G4Decay::G4Decay(), G4DecayWithSpin::G4DecayWithSpin(), G4DNAAttachment::G4DNAAttachment(), G4DNABrownianTransportation::G4DNABrownianTransportation(), G4DNAChargeDecrease::G4DNAChargeDecrease(), G4DNAChargeIncrease::G4DNAChargeIncrease(), G4DNAElastic::G4DNAElastic(), G4DNAElectronSolvatation::G4DNAElectronSolvatation(), G4DNAExcitation::G4DNAExcitation(), G4DNAIonisation::G4DNAIonisation(), G4DNAMolecularDecay::G4DNAMolecularDecay(), G4DNAVibExcitation::G4DNAVibExcitation(), G4eBremsstrahlung::G4eBremsstrahlung(), G4eeToHadrons::G4eeToHadrons(), G4eIonisation::G4eIonisation(), G4eplusAnnihilation::G4eplusAnnihilation(), G4eplusPolarizedAnnihilation::G4eplusPolarizedAnnihilation(), G4ePolarizedIonisation::G4ePolarizedIonisation(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4GammaConversion::G4GammaConversion(), G4GammaConversionToMuons::G4GammaConversionToMuons(), G4HadronicProcess::G4HadronicProcess(), G4hBremsstrahlung::G4hBremsstrahlung(), G4hhIonisation::G4hhIonisation(), G4hIonisation::G4hIonisation(), G4hPairProduction::G4hPairProduction(), G4ionIonisation::G4ionIonisation(), G4ITTransportation::G4ITTransportation(), G4KaonMinusAbsorption::G4KaonMinusAbsorption(), G4KaonMinusAbsorptionAtRest::G4KaonMinusAbsorptionAtRest(), G4mplIonisation::G4mplIonisation(), G4MuBremsstrahlung::G4MuBremsstrahlung(), G4MuElecElastic::G4MuElecElastic(), G4MuElecInelastic::G4MuElecInelastic(), G4MuIonisation::G4MuIonisation(), G4MuonMinusCaptureAtRest::G4MuonMinusCaptureAtRest(), G4MuPairProduction::G4MuPairProduction(), G4NeutronCaptureAtRest::G4NeutronCaptureAtRest(), G4NeutronKiller::G4NeutronKiller(), G4NuclearStopping::G4NuclearStopping(), G4OpAbsorption::G4OpAbsorption(), G4OpBoundaryProcess::G4OpBoundaryProcess(), G4OpMieHG::G4OpMieHG(), G4OpRayleigh::G4OpRayleigh(), G4OpWLS::G4OpWLS(), G4PhotoElectricEffect::G4PhotoElectricEffect(), G4PiMinusAbsorptionAtRest::G4PiMinusAbsorptionAtRest(), G4PionDecayMakeSpin::G4PionDecayMakeSpin(), G4PionMinusAbsorptionAtRest::G4PionMinusAbsorptionAtRest(), G4PionMinusNuclearAtRestChips::G4PionMinusNuclearAtRestChips(), G4PolarizedCompton::G4PolarizedCompton(), G4PolarizedGammaConversion::G4PolarizedGammaConversion(), G4PolarizedPhotoElectricEffect::G4PolarizedPhotoElectricEffect(), G4ProtonAntiProtonAtRestChips::G4ProtonAntiProtonAtRestChips(), G4RadioactiveDecay::G4RadioactiveDecay(), G4RayleighScattering::G4RayleighScattering(), G4Scintillation::G4Scintillation(), G4StepLimiter::G4StepLimiter(), G4SynchrotronRadiation::G4SynchrotronRadiation(), G4SynchrotronRadiationInMat::G4SynchrotronRadiationInMat(), G4TransitionRadiation::G4TransitionRadiation(), G4Transportation::G4Transportation(), G4UnknownDecay::G4UnknownDecay(), G4UserSpecialCuts::G4UserSpecialCuts(), G4VMultipleScattering::G4VMultipleScattering(), and G4Decay::SetExtDecayer().

00404 {
00405  theProcessSubType = value;
00406 }

void G4VProcess::SetProcessType ( G4ProcessType   )  [inline]

Definition at line 391 of file G4VProcess.hh.

References theProcessType.

Referenced by MaxTimeCuts::MaxTimeCuts(), and MinEkineCuts::MinEkineCuts().

00392 {
00393   theProcessType = aType;
00394 }

void G4VProcess::SetVerboseLevel ( G4int  value  )  [inline]

Reimplemented in G4Decay, G4UnknownDecay, G4ITTransportation, G4SynchrotronRadiationInMat, G4RadioactiveDecay, G4CoupledTransportation, and G4Transportation.

Definition at line 408 of file G4VProcess.hh.

References verboseLevel.

Referenced by TLBE< T >::ConstructOp(), G4ProcessTable::DumpInfo(), G4eeToHadrons::G4eeToHadrons(), G4hhIonisation::G4hhIonisation(), G4mplIonisation::G4mplIonisation(), G4VEmProcess::G4VEmProcess(), G4VEnergyLossProcess::G4VEnergyLossProcess(), G4VMultipleScattering::G4VMultipleScattering(), G4CoulombScattering::InitialiseProcess(), G4ProcessTableMessenger::SetNewValue(), and G4ProcessManagerMessenger::SetNewValue().

00409 {
00410   verboseLevel = value;
00411 }

void G4VProcess::StartTracking ( G4Track  )  [virtual]

Reimplemented in G4ImportanceProcess, G4WeightCutOffProcess, G4WeightWindowProcess, G4Decay, G4AdjointProcessEquivalentToDirectProcess, G4ITTransportation, G4VITProcess, G4DNABrownianTransportation, G4DNASecondOrderReaction, G4VEmProcess, G4VEnergyLossProcess, G4VMultipleScattering, G4WrapperProcess, G4FastSimulationManagerProcess, G4ParallelWorldProcess, G4ParallelWorldScoringProcess, G4ScoreSplittingProcess, G4CoupledTransportation, and G4Transportation.

Definition at line 125 of file G4VProcess.cc.

References currentInteractionLength, G4cout, G4endl, theInitialNumberOfInteractionLength, theNumberOfInteractionLengthLeft, theProcessName, and verboseLevel.

Referenced by G4WrapperProcess::StartTracking(), G4Transportation::StartTracking(), G4ITTransportation::StartTracking(), G4DNASecondOrderReaction::StartTracking(), and G4AdjointProcessEquivalentToDirectProcess::StartTracking().

00126 {
00127   currentInteractionLength = -1.0;
00128   theNumberOfInteractionLengthLeft = -1.0;
00129   theInitialNumberOfInteractionLength=-1.0;
00130 #ifdef G4VERBOSE
00131   if (verboseLevel>2) {
00132     G4cout << "G4VProcess::StartTracking() [" << theProcessName << "]" <<G4endl;
00133   }
00134 #endif
00135 }

virtual G4bool G4VProcess::StorePhysicsTable ( const G4ParticleDefinition ,
const G4String ,
G4bool   
) [inline, virtual]

Reimplemented in G4AdjointProcessEquivalentToDirectProcess, G4VEmProcess, G4VEnergyLossProcess, G4VMultipleScattering, and G4WrapperProcess.

Definition at line 231 of file G4VProcess.hh.

Referenced by G4WrapperProcess::StorePhysicsTable(), and G4AdjointProcessEquivalentToDirectProcess::StorePhysicsTable().

00232                                                                {return true;}

void G4VProcess::SubtractNumberOfInteractionLengthLeft ( G4double  previousStepSize  )  [protected]

Reimplemented in G4VITProcess, G4IVContinuousDiscreteProcess, and G4IVRestDiscreteProcess.

Definition at line 98 of file G4VProcess.cc.

References currentInteractionLength, EventMustBeAborted, G4cerr, G4endl, G4Exception(), theNumberOfInteractionLengthLeft, theProcessName, and verboseLevel.

Referenced by G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), and G4Decay::PostStepGetPhysicalInteractionLength().

00100 {
00101   if (currentInteractionLength>0.0) {
00102     theNumberOfInteractionLengthLeft -= previousStepSize/currentInteractionLength;
00103     if(theNumberOfInteractionLengthLeft<0.) {
00104        theNumberOfInteractionLengthLeft=perMillion;
00105     }          
00106 
00107   } else {
00108 #ifdef G4VERBOSE
00109     if (verboseLevel>0) {
00110       G4cerr << "G4VProcess::SubtractNumberOfInteractionLengthLeft()";
00111       G4cerr << " [" << theProcessName << "]" <<G4endl;
00112       G4cerr << " currentInteractionLength = " << currentInteractionLength/cm << " [cm]";
00113       G4cerr << " previousStepSize = " << previousStepSize/cm << " [cm]";
00114       G4cerr << G4endl;
00115     }
00116 #endif
00117     G4String msg = "Negative currentInteractionLength for ";
00118     msg +=      theProcessName;
00119     G4Exception("G4VProcess::SubtractNumberOfInteractionLengthLeft()",
00120                 "ProcMan201",EventMustBeAborted,
00121                 msg);
00122   }
00123 }


Field Documentation

G4ParticleChange G4VProcess::aParticleChange [protected]

Definition at line 289 of file G4VProcess.hh.

Referenced by G4hImpactIonisation::AlongStepDoIt(), G4ErrorEnergyLoss::AlongStepDoIt(), G4ContinuousGainOfEnergy::AlongStepDoIt(), G4QCaptureAtRest::AtRestDoIt(), G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4PiMinusAbsorptionAtRest::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4MuonMinusCaptureAtRest::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4DNAMolecularDecay::DecayIt(), G4HadronicProcess::FillResult(), G4DNAMolecularDecay::G4DNAMolecularDecay(), G4VProcess(), SpecialCuts::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4TransitionRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4StepLimiter::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4QSynchRad::PostStepDoIt(), G4QNGamma::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QIonIonElastic::PostStepDoIt(), G4QInelastic::PostStepDoIt(), G4QElastic::PostStepDoIt(), G4QDiffraction::PostStepDoIt(), G4QCoherentChargeExchange::PostStepDoIt(), G4QAtomicElectronScattering::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpRayleigh::PostStepDoIt(), G4OpMieHG::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), and G4AnnihiToMuPair::PostStepDoIt().

const G4ProcessManager* G4VProcess::aProcessManager [protected]

Definition at line 280 of file G4VProcess.hh.

Referenced by GetProcessManager(), and SetProcessManager().

G4double G4VProcess::currentInteractionLength [protected]

Definition at line 297 of file G4VProcess.hh.

Referenced by G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4ProtonAntiProtonAtRestChips::AtRestGetPhysicalInteractionLength(), G4PionMinusNuclearAtRestChips::AtRestGetPhysicalInteractionLength(), G4PionMinusAbsorptionAtRest::AtRestGetPhysicalInteractionLength(), G4NeutronCaptureAtRest::AtRestGetPhysicalInteractionLength(), G4KaonMinusAbsorption::AtRestGetPhysicalInteractionLength(), G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4AntiProtonAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), G4AntiNeutronAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), EndTracking(), G4Decay::EndTracking(), GetCurrentInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4Decay::PostStepGetPhysicalInteractionLength(), StartTracking(), G4Decay::StartTracking(), and SubtractNumberOfInteractionLengthLeft().

G4bool G4VProcess::enableAlongStepDoIt [protected]

Definition at line 351 of file G4VProcess.hh.

Referenced by G4DNAMolecularDecay::G4DNAMolecularDecay(), G4ITTransportation::G4ITTransportation(), G4IVRestDiscreteProcess::G4IVRestDiscreteProcess(), G4NuclearStopping::G4NuclearStopping(), G4VDiscreteProcess::G4VDiscreteProcess(), G4VITRestDiscreteProcess::G4VITRestDiscreteProcess(), G4VITRestProcess::G4VITRestProcess(), G4VRestDiscreteProcess::G4VRestDiscreteProcess(), G4VRestProcess::G4VRestProcess(), and isAlongStepDoItIsEnabled().

G4bool G4VProcess::enableAtRestDoIt [protected]

Definition at line 350 of file G4VProcess.hh.

Referenced by G4DNAMolecularDecay::G4DNAMolecularDecay(), G4eplusAnnihilation::G4eplusAnnihilation(), G4eplusPolarizedAnnihilation::G4eplusPolarizedAnnihilation(), G4HadronStoppingProcess::G4HadronStoppingProcess(), G4ITTransportation::G4ITTransportation(), G4IVContinuousDiscreteProcess::G4IVContinuousDiscreteProcess(), G4VContinuousDiscreteProcess::G4VContinuousDiscreteProcess(), G4VContinuousProcess::G4VContinuousProcess(), G4VDiscreteProcess::G4VDiscreteProcess(), and isAtRestDoItIsEnabled().

G4bool G4VProcess::enablePostStepDoIt [protected]

Definition at line 352 of file G4VProcess.hh.

Referenced by G4DNAMolecularDecay::G4DNAMolecularDecay(), G4HadronStoppingProcess::G4HadronStoppingProcess(), G4ITTransportation::G4ITTransportation(), G4NuclearStopping::G4NuclearStopping(), G4VContinuousProcess::G4VContinuousProcess(), G4VITRestProcess::G4VITRestProcess(), G4VRestContinuousProcess::G4VRestContinuousProcess(), G4VRestProcess::G4VRestProcess(), and isPostStepDoItIsEnabled().

G4VParticleChange* G4VProcess::pParticleChange [protected]

Definition at line 283 of file G4VProcess.hh.

Referenced by G4VMultipleScattering::AddEmModel(), G4VEnergyLossProcess::AddEmModel(), G4VEmProcess::AddEmModel(), G4WeightWindowProcess::AlongStepDoIt(), G4WeightCutOffProcess::AlongStepDoIt(), G4VRestContinuousProcess::AlongStepDoIt(), G4VRestContinuousDiscreteProcess::AlongStepDoIt(), G4VContinuousProcess::AlongStepDoIt(), G4VContinuousDiscreteProcess::AlongStepDoIt(), G4ParallelWorldScoringProcess::AlongStepDoIt(), G4ParallelWorldProcess::AlongStepDoIt(), G4IVContinuousDiscreteProcess::AlongStepDoIt(), G4ImportanceProcess::AlongStepDoIt(), G4VRestProcess::AtRestDoIt(), G4VRestDiscreteProcess::AtRestDoIt(), G4VRestContinuousProcess::AtRestDoIt(), G4VRestContinuousDiscreteProcess::AtRestDoIt(), G4VITRestProcess::AtRestDoIt(), G4VITRestDiscreteProcess::AtRestDoIt(), G4ScoreSplittingProcess::AtRestDoIt(), G4ParallelWorldScoringProcess::AtRestDoIt(), G4ParallelWorldProcess::AtRestDoIt(), G4IVRestDiscreteProcess::AtRestDoIt(), G4Decay::G4Decay(), G4DNAMolecularDecay::G4DNAMolecularDecay(), G4ImportanceProcess::G4ImportanceProcess(), G4ITTransportation::G4ITTransportation(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess(), G4RadioactiveDecay::G4RadioactiveDecay(), G4ScoreSplittingProcess::G4ScoreSplittingProcess(), G4UnknownDecay::G4UnknownDecay(), G4VEmProcess::G4VEmProcess(), G4VEnergyLossProcess::G4VEnergyLossProcess(), G4VMultipleScattering::G4VMultipleScattering(), G4VProcess(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4WeightCutOffProcess::G4WeightCutOffProcess(), G4WeightWindowProcess::G4WeightWindowProcess(), G4VTransitionRadiation::PostStepDoIt(), G4VRestDiscreteProcess::PostStepDoIt(), G4VRestContinuousDiscreteProcess::PostStepDoIt(), G4VITRestDiscreteProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VDiscreteProcess::PostStepDoIt(), G4VContinuousDiscreteProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4ParallelWorldScoringProcess::PostStepDoIt(), G4ParallelWorldProcess::PostStepDoIt(), G4NeutronKiller::PostStepDoIt(), G4IVRestDiscreteProcess::PostStepDoIt(), G4IVContinuousDiscreteProcess::PostStepDoIt(), and G4Cerenkov::PostStepDoIt().

G4double G4VProcess::theInitialNumberOfInteractionLength [protected]

Definition at line 300 of file G4VProcess.hh.

Referenced by ClearNumberOfInteractionLengthLeft(), EndTracking(), GetTotalNumberOfInteractionLengthTraversed(), ResetNumberOfInteractionLengthLeft(), and StartTracking().

G4double G4VProcess::theNumberOfInteractionLengthLeft [protected]

Definition at line 293 of file G4VProcess.hh.

Referenced by G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4ProtonAntiProtonAtRestChips::AtRestGetPhysicalInteractionLength(), G4PionMinusNuclearAtRestChips::AtRestGetPhysicalInteractionLength(), G4PionMinusAbsorptionAtRest::AtRestGetPhysicalInteractionLength(), G4NeutronCaptureAtRest::AtRestGetPhysicalInteractionLength(), G4KaonMinusAbsorption::AtRestGetPhysicalInteractionLength(), G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4Decay::AtRestGetPhysicalInteractionLength(), G4AntiProtonAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), G4AntiNeutronAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), ClearNumberOfInteractionLengthLeft(), EndTracking(), GetNumberOfInteractionLengthLeft(), GetTotalNumberOfInteractionLengthTraversed(), G4VEnergyLossProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4Decay::PostStepGetPhysicalInteractionLength(), ResetNumberOfInteractionLengthLeft(), StartTracking(), G4VEnergyLossProcess::StartTracking(), G4VEmProcess::StartTracking(), and SubtractNumberOfInteractionLengthLeft().

G4String G4VProcess::thePhysicsTableFileName [protected]

Definition at line 338 of file G4VProcess.hh.

Referenced by GetPhysicsTableFileName().

G4double G4VProcess::thePILfactor [protected]

Definition at line 346 of file G4VProcess.hh.

Referenced by AtRestGPIL(), GetPILfactor(), PostStepGPIL(), and SetPILfactor().

G4String G4VProcess::theProcessName [protected]

Definition at line 335 of file G4VProcess.hh.

Referenced by DumpInfo(), EndTracking(), G4WeightCutOffProcess::GetName(), G4ImportanceProcess::GetName(), GetPhysicsTableFileName(), GetProcessName(), G4WrapperProcess::RegisterProcess(), StartTracking(), SubtractNumberOfInteractionLengthLeft(), and G4VITProcess::SubtractNumberOfInteractionLengthLeft().

G4int G4VProcess::theProcessSubType [protected]

Definition at line 343 of file G4VProcess.hh.

Referenced by DumpInfo(), GetProcessSubType(), and SetProcessSubType().

G4ProcessType G4VProcess::theProcessType [protected]

Definition at line 340 of file G4VProcess.hh.

Referenced by DumpInfo(), G4AdjointProcessEquivalentToDirectProcess::G4AdjointProcessEquivalentToDirectProcess(), GetProcessType(), G4WrapperProcess::RegisterProcess(), and SetProcessType().

G4int G4VProcess::verboseLevel [protected]

Reimplemented in G4Decay.

Definition at line 368 of file G4VProcess.hh.

Referenced by G4VEnergyLossProcess::ActivateForcedInteraction(), G4VEmProcess::ActivateForcedInteraction(), G4VEnergyLossProcess::ActivateSecondaryBiasing(), G4VEmProcess::ActivateSecondaryBiasing(), G4VEnergyLossProcess::AddCollaborativeProcess(), G4VEnergyLossProcess::AlongStepDoIt(), G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4IVContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4PiMinusAbsorptionAtRest::AtRestDoIt(), G4ParallelWorldScoringProcess::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestProcess::AtRestGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4ProtonAntiProtonAtRestChips::AtRestGetPhysicalInteractionLength(), G4PionMinusNuclearAtRestChips::AtRestGetPhysicalInteractionLength(), G4PionMinusAbsorptionAtRest::AtRestGetPhysicalInteractionLength(), G4NeutronCaptureAtRest::AtRestGetPhysicalInteractionLength(), G4KaonMinusAbsorption::AtRestGetPhysicalInteractionLength(), G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4AntiProtonAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), G4AntiNeutronAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), G4VXTRenergyLoss::BuildAngleForEnergyBank(), G4VXTRenergyLoss::BuildAngleTable(), G4VEnergyLossProcess::BuildDEDXTable(), G4VXTRenergyLoss::BuildEnergyTable(), G4VXTRenergyLoss::BuildGlobalAngleTable(), G4VEnergyLossProcess::BuildLambdaTable(), G4VXTRenergyLoss::BuildPhysicsTable(), G4VMultipleScattering::BuildPhysicsTable(), G4VEnergyLossProcess::BuildPhysicsTable(), G4VEmProcess::BuildPhysicsTable(), G4SynchrotronRadiation::BuildPhysicsTable(), G4hImpactIonisation::BuildPhysicsTable(), G4DNABrownianTransportation::BuildPhysicsTable(), G4ChargeExchangeProcess::BuildPhysicsTable(), EndTracking(), G4AntiNeutronAnnihilationAtRest::G4AntiNeutronAnnihilationAtRest(), G4AntiProtonAnnihilationAtRest::G4AntiProtonAnnihilationAtRest(), G4Cerenkov::G4Cerenkov(), G4ChargeExchangeProcess::G4ChargeExchangeProcess(), G4DNABrownianTransportation::G4DNABrownianTransportation(), G4DNAMolecularDecay::G4DNAMolecularDecay(), G4ePolarizedIonisation::G4ePolarizedIonisation(), G4ErrorEnergyLoss::G4ErrorEnergyLoss(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4ImportanceProcess::G4ImportanceProcess(), G4KaonMinusAbsorption::G4KaonMinusAbsorption(), G4KaonMinusAbsorptionAtRest::G4KaonMinusAbsorptionAtRest(), G4NeutronCaptureAtRest::G4NeutronCaptureAtRest(), G4OpAbsorption::G4OpAbsorption(), G4OpBoundaryProcess::G4OpBoundaryProcess(), G4OpMieHG::G4OpMieHG(), G4OpRayleigh::G4OpRayleigh(), G4OpWLS::G4OpWLS(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess(), G4PiMinusAbsorptionAtRest::G4PiMinusAbsorptionAtRest(), G4PionMinusAbsorptionAtRest::G4PionMinusAbsorptionAtRest(), G4QAtomicElectronScattering::G4QAtomicElectronScattering(), G4QCaptureAtRest::G4QCaptureAtRest(), G4QCoherentChargeExchange::G4QCoherentChargeExchange(), G4QDiffraction::G4QDiffraction(), G4QDiscProcessMixer::G4QDiscProcessMixer(), G4QElastic::G4QElastic(), G4QInelastic::G4QInelastic(), G4QIonIonElastic::G4QIonIonElastic(), G4QLowEnergy::G4QLowEnergy(), G4QNGamma::G4QNGamma(), G4Scintillation::G4Scintillation(), G4ScoreSplittingProcess::G4ScoreSplittingProcess(), G4StepLimiter::G4StepLimiter(), G4StrawTubeXTRadiator::G4StrawTubeXTRadiator(), G4SynchrotronRadiation::G4SynchrotronRadiation(), G4TransparentRegXTRadiator::G4TransparentRegXTRadiator(), G4UserSpecialCuts::G4UserSpecialCuts(), G4VLowEnergyDiscretePhotonProcess::G4VLowEnergyDiscretePhotonProcess(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4WeightCutOffProcess::G4WeightCutOffProcess(), G4WeightWindowProcess::G4WeightWindowProcess(), G4VXTRenergyLoss::GetAngleVector(), G4ChargeExchangeProcess::GetElementCrossSection(), G4VXTRenergyLoss::GetGasZmuProduct(), G4VXTRenergyLoss::GetMeanFreePath(), G4SynchrotronRadiation::GetMeanFreePath(), G4PolarizedCompton::GetMeanFreePath(), G4eplusPolarizedAnnihilation::GetMeanFreePath(), G4VXTRenergyLoss::GetNumberOfPhotons(), G4VXTRenergyLoss::GetPlateZmuProduct(), G4SynchrotronRadiation::GetRandomEnergySR(), GetVerboseLevel(), G4hhIonisation::InitialiseEnergyLossProcess(), G4eeToHadrons::InitialiseProcess(), MaxTimeCuts::MaxTimeCuts(), MinEkineCuts::MinEkineCuts(), G4WHadronElasticProcess::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4ParallelWorldScoringProcess::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpRayleigh::PostStepDoIt(), G4OpMieHG::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4DNASecondOrderReaction::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4PolarizedCompton::PostStepGetPhysicalInteractionLength(), G4eplusPolarizedAnnihilation::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4VMultipleScattering::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), G4VEmProcess::PreparePhysicsTable(), G4VMultipleScattering::PrintInfoDefinition(), G4VEnergyLossProcess::PrintInfoDefinition(), G4VEmProcess::PrintInfoDefinition(), G4VEnergyLossProcess::RetrievePhysicsTable(), G4VEmProcess::RetrievePhysicsTable(), G4VEnergyLossProcess::SetCrossSectionBiasingFactor(), G4VEmProcess::SetCrossSectionBiasingFactor(), G4VEnergyLossProcess::SetInverseRangeTable(), G4NeutronKiller::SetKinEnergyLimit(), G4VEnergyLossProcess::SetLambdaTable(), G4VEnergyLossProcess::SetRangeTableForLoss(), G4VEnergyLossProcess::SetSecondaryRangeTable(), G4VEnergyLossProcess::SetSubLambdaTable(), G4NeutronKiller::SetTimeLimit(), SetVerboseLevel(), G4FastSimulationManagerProcess::SetWorldVolume(), SpecialCuts::SpecialCuts(), G4XTRRegularRadModel::SpectralXTRdEdx(), G4TransparentRegXTRadiator::SpectralXTRdEdx(), G4RegularXTRadiator::SpectralXTRdEdx(), StartTracking(), G4VMultipleScattering::StorePhysicsTable(), G4VEnergyLossProcess::StorePhysicsTable(), SubtractNumberOfInteractionLengthLeft(), G4VITProcess::SubtractNumberOfInteractionLengthLeft(), G4VEmProcess::~G4VEmProcess(), G4VEnergyLossProcess::~G4VEnergyLossProcess(), and G4VMultipleScattering::~G4VMultipleScattering().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:53 2013 for Geant4 by  doxygen 1.4.7