Geant4-11
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Attributes
G4DNAMolecularDissociation Class Reference

#include <G4DNAMolecularDissociation.hh>

Inheritance diagram for G4DNAMolecularDissociation:
G4VITRestDiscreteProcess G4VITProcess G4VProcess

Public Types

using Displacer = G4VMolecularDissociationDisplacer
 
using Species = const G4MoleculeDefinition
 

Public Member Functions

G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &step) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition) override
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void DumpInfo () const
 
virtual void EndTracking ()
 
 G4DNAMolecularDissociation ()=delete
 
 G4DNAMolecularDissociation (const G4DNAMolecularDissociation &right)=delete
 
 G4DNAMolecularDissociation (const G4String &processName="DNAMolecularDecay", G4ProcessType type=fDecay)
 
G4double GetCurrentInteractionLength () const
 
DisplacerGetDisplacer (Species *)
 
G4double GetInteractionTimeLeft ()
 
const G4VProcessGetMasterProcess () const
 
G4double GetNumberOfInteractionLengthLeft () const
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
G4double GetPILfactor () const
 
size_t GetProcessID () const
 
virtual const G4ProcessManagerGetProcessManager ()
 
const G4StringGetProcessName () const
 
G4shared_ptr< G4ProcessState_LockGetProcessState ()
 
G4int GetProcessSubType () const
 
G4ProcessType GetProcessType () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4int GetVerboseLevel () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
G4bool operator!= (const G4VITProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4DNAMolecularDissociationoperator= (const G4DNAMolecularDissociation &right)=delete
 
G4bool operator== (const G4VITProcess &right) const
 
G4bool operator== (const G4VProcess &right) const
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &step) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 
virtual void ProcessDescription (std::ostream &outfile) const
 
G4bool ProposesTimeStep () const
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
void ResetProcessState ()
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
void SetDisplacer (Species *, Displacer *)
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
void SetPILfactor (G4double value)
 
virtual void SetProcessManager (const G4ProcessManager *)
 
void SetProcessState (G4shared_ptr< G4ProcessState_Lock > aProcInfo)
 
void SetProcessSubType (G4int)
 
void SetProcessType (G4ProcessType)
 
void SetVerbose (G4int)
 
void SetVerboseLevel (G4int value)
 
virtual void StartTracking (G4Track *)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
 ~G4DNAMolecularDissociation () override
 

Static Public Member Functions

static const size_t & GetMaxProcessIndex ()
 
static const G4StringGetProcessTypeName (G4ProcessType)
 

Protected Member Functions

virtual void ClearInteractionTimeLeft ()
 
virtual void ClearNumberOfInteractionLengthLeft ()
 
void CreateInfo ()
 
virtual G4VParticleChangeDecayIt (const G4Track &, const G4Step &)
 
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *) override
 
G4double GetMeanLifeTime (const G4Track &, G4ForceCondition *) override
 
template<typename T >
T * GetState ()
 
G4bool InstantiateProcessState ()
 
void RetrieveProcessInfo ()
 
void SetInstantiateProcessState (G4bool flag)
 
virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 

Protected Attributes

G4ParticleChange aParticleChange
 
const G4ProcessManageraProcessManager = nullptr
 
G4bool enableAlongStepDoIt = true
 
G4bool enableAtRestDoIt = true
 
G4bool enablePostStepDoIt = true
 
G4bool fProposesTimeStep
 
G4shared_ptr< G4ProcessStatefpState
 
G4VParticleChangepParticleChange = nullptr
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String thePhysicsTableFileName
 
G4double thePILfactor = 1.0
 
G4String theProcessName
 
G4int theProcessSubType = -1
 
G4ProcessType theProcessType = fNotDefined
 
G4int verboseLevel = 0
 

Private Types

using DisplacementMap = std::map< Species *, std::unique_ptr< Displacer > >
 

Private Attributes

G4doublecurrentInteractionLength
 
G4bool fDecayAtFixedTime
 
DisplacementMap fDisplacementMap
 
G4bool fInstantiateProcessState
 
static size_t * fNbProcess = 0
 
size_t fProcessID
 
G4ProcessTablefProcessTable = nullptr
 
G4int fVerbose
 
G4VProcessmasterProcessShadow = nullptr
 
G4doubletheInteractionTimeLeft
 
G4doubletheNumberOfInteractionLengthLeft
 

Detailed Description

G4DNAMolecularDissociation should be called only for molecules. It will dissociate the molecules using the decay associated to this molecule and if a displacement scheme has been registered, it will place the products to the expected position.

Definition at line 63 of file G4DNAMolecularDissociation.hh.

Member Typedef Documentation

◆ DisplacementMap

using G4DNAMolecularDissociation::DisplacementMap = std::map<Species*, std::unique_ptr<Displacer> >
private

Definition at line 104 of file G4DNAMolecularDissociation.hh.

◆ Displacer

Definition at line 74 of file G4DNAMolecularDissociation.hh.

◆ Species

Definition at line 73 of file G4DNAMolecularDissociation.hh.

Constructor & Destructor Documentation

◆ G4DNAMolecularDissociation() [1/3]

G4DNAMolecularDissociation::G4DNAMolecularDissociation ( const G4String processName = "DNAMolecularDecay",
G4ProcessType  type = fDecay 
)

Definition at line 50 of file G4DNAMolecularDissociation.cc.

53 : G4VITRestDiscreteProcess(processName, type)
54{
55 // set Process Sub Type
57 enableAlongStepDoIt = false;
58 enablePostStepDoIt = true;
59 enableAtRestDoIt = true;
60
61 fVerbose = 0;
62
63#ifdef G4VERBOSE
64 if (verboseLevel > 1)
65 {
66 G4cout << "G4MolecularDissociationProcess constructor " << " Name:"
67 << processName << G4endl;
68 }
69#endif
70
72
73 fDecayAtFixedTime = true;
74 fProposesTimeStep = true;
75}
@ fLowEnergyMolecularDecay
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4bool fProposesTimeStep
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:359
G4int verboseLevel
Definition: G4VProcess.hh:356
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:361
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321

References G4VProcess::aParticleChange, G4VProcess::enableAlongStepDoIt, G4VProcess::enableAtRestDoIt, G4VProcess::enablePostStepDoIt, fDecayAtFixedTime, fLowEnergyMolecularDecay, G4VITProcess::fProposesTimeStep, fVerbose, G4cout, G4endl, G4VProcess::pParticleChange, G4VProcess::SetProcessSubType(), and G4VProcess::verboseLevel.

◆ G4DNAMolecularDissociation() [2/3]

G4DNAMolecularDissociation::G4DNAMolecularDissociation ( )
delete

◆ G4DNAMolecularDissociation() [3/3]

G4DNAMolecularDissociation::G4DNAMolecularDissociation ( const G4DNAMolecularDissociation right)
delete

◆ ~G4DNAMolecularDissociation()

G4DNAMolecularDissociation::~G4DNAMolecularDissociation ( )
overridedefault

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4VITRestDiscreteProcess::AlongStepDoIt ( const G4Track ,
const G4Step  
)
inlineoverridevirtualinherited

Implements G4VProcess.

Definition at line 92 of file G4VITRestDiscreteProcess.hh.

93 {
94 return nullptr;
95 }

◆ AlongStepGetPhysicalInteractionLength()

G4double G4VITRestDiscreteProcess::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4double  ,
G4double ,
G4GPILSelection  
)
inlineoverridevirtualinherited

Implements G4VProcess.

Definition at line 82 of file G4VITRestDiscreteProcess.hh.

87 {
88 return -1.0;
89 }

◆ AlongStepGPIL()

G4double G4VProcess::AlongStepGPIL ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
inlineinherited

Definition at line 461 of file G4VProcess.hh.

466{
467 return AlongStepGetPhysicalInteractionLength(track, previousStepSize,
468 currentMinimumStep, proposedSafety, selection);
469}
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0

References G4VProcess::AlongStepGetPhysicalInteractionLength().

Referenced by G4SteppingManager::DefinePhysicalStepLength(), and G4ITStepProcessor::DoDefinePhysicalStepLength().

◆ AtRestDoIt()

G4VParticleChange * G4DNAMolecularDissociation::AtRestDoIt ( const G4Track track,
const G4Step step 
)
overridevirtual

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 315 of file G4DNAMolecularDissociation.cc.

317{
320 return DecayIt(track, step);
321}
virtual G4VParticleChange * DecayIt(const G4Track &, const G4Step &)
virtual void ClearInteractionTimeLeft()
virtual void ClearNumberOfInteractionLengthLeft()

References G4VITProcess::ClearInteractionTimeLeft(), G4VITProcess::ClearNumberOfInteractionLengthLeft(), and DecayIt().

Referenced by PostStepDoIt().

◆ AtRestGetPhysicalInteractionLength()

G4double G4DNAMolecularDissociation::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
overridevirtual

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 325 of file G4DNAMolecularDissociation.cc.

327{
329 {
330 return GetMeanLifeTime(track, condition);
331 }
332
334}
G4double condition(const G4ErrorSymMatrix &m)
G4double GetMeanLifeTime(const G4Track &, G4ForceCondition *) override
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override

References G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), condition(), fDecayAtFixedTime, and GetMeanLifeTime().

◆ AtRestGPIL()

G4double G4VProcess::AtRestGPIL ( const G4Track track,
G4ForceCondition condition 
)
inlineinherited

Definition at line 472 of file G4VProcess.hh.

474{
476}
G4double thePILfactor
Definition: G4VProcess.hh:352
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)=0

References G4VProcess::AtRestGetPhysicalInteractionLength(), condition(), and G4VProcess::thePILfactor.

Referenced by G4ITStepProcessor::GetAtRestIL(), and G4SteppingManager::InvokeAtRestDoItProcs().

◆ BuildPhysicsTable()

virtual void G4VITProcess::BuildPhysicsTable ( const G4ParticleDefinition )
inlinevirtualinherited

◆ BuildWorkerPhysicsTable()

void G4VProcess::BuildWorkerPhysicsTable ( const G4ParticleDefinition part)
virtualinherited

Reimplemented in G4BiasingProcessInterface.

Definition at line 200 of file G4VProcess.cc.

201{
202 BuildPhysicsTable(part);
203}
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:187

References G4VProcess::BuildPhysicsTable().

Referenced by G4BiasingProcessInterface::BuildWorkerPhysicsTable().

◆ ClearInteractionTimeLeft()

void G4VITProcess::ClearInteractionTimeLeft ( )
inlineprotectedvirtualinherited

Definition at line 258 of file G4VITProcess.hh.

259{
260 fpState->theInteractionTimeLeft = -1.0;
261}
G4shared_ptr< G4ProcessState > fpState

References G4VITProcess::fpState.

Referenced by G4VITRestProcess::AtRestDoIt(), G4DNAElectronHoleRecombination::AtRestDoIt(), and AtRestDoIt().

◆ ClearNumberOfInteractionLengthLeft()

void G4VITProcess::ClearNumberOfInteractionLengthLeft ( )
inlineprotectedvirtualinherited

◆ CreateInfo()

void G4VITProcess::CreateInfo ( )
protectedinherited

◆ DecayIt()

G4VParticleChange * G4DNAMolecularDissociation::DecayIt ( const G4Track track,
const G4Step  
)
protectedvirtual

Definition at line 116 of file G4DNAMolecularDissociation.cc.

118{
120 auto pMotherMolecule = GetMolecule(track);
121 auto pMotherMoleculeDefinition = pMotherMolecule->GetDefinition();
122
123 if (pMotherMoleculeDefinition->GetDecayTable())
124 {
125 const auto pDissociationChannels = pMotherMolecule->GetDissociationChannels();
126
127 if (pDissociationChannels == nullptr)
128 {
129 G4ExceptionDescription exceptionDescription;
130 pMotherMolecule->PrintState();
131 exceptionDescription << "No decay channel was found for the molecule : "
132 << pMotherMolecule->GetName() << G4endl;
133 G4Exception("G4DNAMolecularDissociation::DecayIt",
134 "G4DNAMolecularDissociation::NoDecayChannel",
136 exceptionDescription);
137 return &aParticleChange;
138 }
139
140 auto decayVectorSize = pDissociationChannels->size();
141 G4double RdmValue = G4UniformRand();
142
143 const G4MolecularDissociationChannel* pDecayChannel = nullptr;
144 size_t i = 0;
145 do
146 {
147 pDecayChannel = (*pDissociationChannels)[i];
148 if (RdmValue < pDecayChannel->GetProbability())
149 {
150 break;
151 }
152 RdmValue -= pDecayChannel->GetProbability();
153 i++;
154 } while (i < decayVectorSize);
155
156 G4double decayEnergy = pDecayChannel->GetEnergy();
157 auto nbProducts = pDecayChannel->GetNbProducts();
158
159 if (decayEnergy > 0.)
160 {
162 }
163
164 if (nbProducts)
165 {
166 std::vector<G4ThreeVector> productsDisplacement(nbProducts);
167 G4ThreeVector motherMoleculeDisplacement;
168
169 auto it = fDisplacementMap.find(pMotherMoleculeDefinition);
170
171 if (it != fDisplacementMap.end())
172 {
173 auto pDisplacer = it->second.get();
174 productsDisplacement = pDisplacer->GetProductsDisplacement(pDecayChannel);
175 motherMoleculeDisplacement =
176 pDisplacer->GetMotherMoleculeDisplacement(pDecayChannel);
177 }
178 else
179 {
181 errMsg << "No G4MolecularDecayProcess::theDecayDisplacementMap["
182 << pMotherMolecule->GetName() + "]";
183 G4Exception("G4MolecularDecayProcess::DecayIt",
184 "DNAMolecularDecay001",
186 errMsg);
187 }
188
190
191#ifdef G4VERBOSE
192 if (fVerbose)
193 {
194 G4cout << "Decay Process : " << pMotherMolecule->GetName()
195 << " (trackID :" << track.GetTrackID() << ") "
196 << pDecayChannel->GetName() << G4endl;
197 }
198#endif
199
201
202 for (G4int j = 0; j < nbProducts; j++)
203 {
204 auto pProduct = new G4Molecule(pDecayChannel->GetProduct(j));
205
206 G4ThreeVector displacement = motherMoleculeDisplacement + productsDisplacement[j];
207 double mag_displacement = displacement.mag();
208 G4ThreeVector displacement_direction = displacement / (mag_displacement + 1e-30);
209
210 double prNewSafety = DBL_MAX;
211
212 //double step =
213 pNavigator->CheckNextStep(track.GetPosition(),
214 displacement_direction,
215 mag_displacement,
216 prNewSafety);
217
218 //if(prNewSafety < mag_displacement || step < mag_displacement)
219 mag_displacement = std::min(prNewSafety * 0.8, mag_displacement);
220
221 G4ThreeVector product_pos = track.GetPosition()
222 + displacement_direction * mag_displacement;
223
224 const G4AffineTransform& transform = pNavigator->GetGlobalToLocalTransform();
225
226 G4ThreeVector localPoint = transform.TransformPoint(product_pos); //track.GetPosition());
227
228 if (track.GetTouchable()->GetSolid()->Inside(localPoint) !=
230 {
232 ED << "The decayed product is outside of the volume : "
233 << track.GetTouchable()->GetVolume()->GetName() << G4endl;
234 G4Exception("G4DNAMolecularDissociation::DecayIt()",
235 "OUTSIDE_OF_MOTHER_VOLUME",
236 JustWarning, ED);
237 }
238
239 auto pSecondary = pProduct->BuildTrack(track.GetGlobalTime(), product_pos);
240
241 pSecondary->SetTrackStatus(fAlive);
242#ifdef G4VERBOSE
243 if (fVerbose)
244 {
245 G4cout << "Product : " << pProduct->GetName() << G4endl;
246 }
247#endif
248 // add the secondary track in the List
249 aParticleChange.G4VParticleChange::AddSecondary(pSecondary);
250 }
251#ifdef G4VERBOSE
252 if (fVerbose)
253 {
254 G4cout << "-------------" << G4endl;
255 }
256#endif
257 }
258 //DEBUG
259 else if (fVerbose && decayEnergy)
260 {
261 G4cout << "No products for this channel" << G4endl;
262 G4cout << "-------------" << G4endl;
263 }
264 /*
265 else if(!decayEnergy && !nbProducts)
266 {
267 G4ExceptionDescription errMsg;
268 errMsg << "There is no products and no energy specified in the molecular "
269 "dissociation channel";
270 G4Exception("G4MolecularDissociationProcess::DecayIt",
271 "DNAMolecularDissociation002",
272 FatalErrorInArgument,
273 errMsg);
274 }
275 */
276 }
277
279
280 return &aParticleChange;
281}
@ JustWarning
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
@ fAlive
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
double mag() const
static G4ITTransportationManager * GetTransportationManager()
G4ITNavigator * GetNavigatorForTracking() const
virtual void Initialize(const G4Track &)
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
const G4VTouchable * GetTouchable() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual G4VSolid * GetSolid(G4int depth=0) const
Definition: G4VTouchable.cc:42
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34
@ kInside
Definition: geomdefs.hh:70
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4bool transform(G4String &input, const G4String &type)
#define DBL_MAX
Definition: templates.hh:62

References G4VProcess::aParticleChange, DBL_MAX, fAlive, FatalErrorInArgument, FatalException, fDisplacementMap, fStopAndKill, fVerbose, G4cout, G4endl, G4Exception(), G4UniformRand, G4MolecularDissociationChannel::GetEnergy(), G4Track::GetGlobalTime(), GetMolecule(), G4VPhysicalVolume::GetName(), G4MolecularDissociationChannel::GetName(), G4ITTransportationManager::GetNavigatorForTracking(), G4MolecularDissociationChannel::GetNbProducts(), G4Track::GetPosition(), G4MolecularDissociationChannel::GetProbability(), G4MolecularDissociationChannel::GetProduct(), G4VTouchable::GetSolid(), G4Track::GetTouchable(), G4Track::GetTrackID(), G4ITTransportationManager::GetTransportationManager(), G4VTouchable::GetVolume(), G4ParticleChange::Initialize(), G4VSolid::Inside(), JustWarning, kInside, CLHEP::Hep3Vector::mag(), G4INCL::Math::min(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4VParticleChange::SetNumberOfSecondaries(), and G4coutFormatters::anonymous_namespace{G4coutFormatters.cc}::transform().

Referenced by AtRestDoIt().

◆ DumpInfo()

void G4VProcess::DumpInfo ( ) const
virtualinherited

Reimplemented in G4AdjointAlongStepWeightCorrection, G4AdjointForcedInteractionForGamma, G4AdjointhMultipleScattering, G4ContinuousGainOfEnergy, G4eAdjointMultipleScattering, G4eInverseBremsstrahlung, G4eInverseCompton, G4eInverseIonisation, G4InversePEEffect, G4IonInverseIonisation, G4PolarizedAnnihilation, G4PolarizedBremsstrahlung, G4PolarizedCompton, G4PolarizedGammaConversion, G4PolarizedIonisation, G4PolarizedPhotoElectric, G4Cerenkov, G4ForwardXrayTR, G4GammaXTRadiator, G4GaussXTRadiator, G4RegularXTRadiator, G4Scintillation, G4StrawTubeXTRadiator, G4SynchrotronRadiation, G4TransitionRadiation, G4TransparentRegXTRadiator, G4VTransitionRadiation, G4VXTRenergyLoss, G4XTRGammaRadModel, G4XTRRegularRadModel, and G4XTRTransparentRegRadModel.

Definition at line 167 of file G4VProcess.cc.

168{
169 G4cout << "Process Name " << theProcessName ;
170 G4cout << " : Type[" << GetProcessTypeName(theProcessType) << "]";
171 G4cout << " : SubType[" << theProcessSubType << "]"<< G4endl;
172}
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:134
G4ProcessType theProcessType
Definition: G4VProcess.hh:346
G4int theProcessSubType
Definition: G4VProcess.hh:349
G4String theProcessName
Definition: G4VProcess.hh:341

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

Referenced by G4ProcessTable::DumpInfo(), export_G4VProcess(), G4Scintillation::ProcessDescription(), G4Cerenkov::ProcessDescription(), and G4ProcessManagerMessenger::SetNewValue().

◆ EndTracking()

void G4VProcess::EndTracking ( )
virtualinherited

◆ GetCurrentInteractionLength()

G4double G4VProcess::GetCurrentInteractionLength ( ) const
inlineinherited

◆ GetDisplacer()

G4DNAMolecularDissociation::Displacer * G4DNAMolecularDissociation::GetDisplacer ( Species pSpecies)

Definition at line 292 of file G4DNAMolecularDissociation.cc.

293{
294 return fDisplacementMap[pSpecies].get();
295}

References fDisplacementMap.

◆ GetInteractionTimeLeft()

G4double G4VITProcess::GetInteractionTimeLeft ( )
inlineinherited

Definition at line 273 of file G4VITProcess.hh.

274{
275 if (fpState) return fpState->theInteractionTimeLeft;
276
277 return -1;
278}

References G4VITProcess::fpState.

Referenced by G4ITStepProcessor::DoDefinePhysicalStepLength().

◆ GetMasterProcess()

const G4VProcess * G4VProcess::GetMasterProcess ( ) const
inlineinherited

◆ GetMaxProcessIndex()

const size_t & G4VITProcess::GetMaxProcessIndex ( )
inlinestaticinherited

Definition at line 285 of file G4VITProcess.hh.

286{
287 if (!fNbProcess) fNbProcess = new size_t(0);
288 return *fNbProcess;
289}
static size_t * fNbProcess

References G4VITProcess::fNbProcess.

Referenced by G4TrackingInformation::GetProcessState().

◆ GetMeanFreePath()

G4double G4DNAMolecularDissociation::GetMeanFreePath ( const G4Track ,
G4double  ,
G4ForceCondition  
)
overrideprotectedvirtual

Implements G4VITRestDiscreteProcess.

Definition at line 346 of file G4DNAMolecularDissociation.cc.

349{
350 return 0;
351}

◆ GetMeanLifeTime()

G4double G4DNAMolecularDissociation::GetMeanLifeTime ( const G4Track track,
G4ForceCondition  
)
overrideprotectedvirtual

Implements G4VITRestDiscreteProcess.

Definition at line 107 of file G4DNAMolecularDissociation.cc.

109{
110 G4double output = GetMolecule(track)->GetDecayTime() - track.GetProperTime();
111 return output > 0. ? output : 0.;
112}
G4double GetDecayTime() const
Definition: G4Molecule.cc:474
G4double GetProperTime() const

References G4Molecule::GetDecayTime(), GetMolecule(), and G4Track::GetProperTime().

Referenced by AtRestGetPhysicalInteractionLength().

◆ GetNumberOfInteractionLengthLeft()

G4double G4VProcess::GetNumberOfInteractionLengthLeft ( ) const
inlineinherited

Definition at line 431 of file G4VProcess.hh.

432{
434}

References G4VProcess::theNumberOfInteractionLengthLeft.

◆ GetPhysicsTableFileName()

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

Definition at line 181 of file G4VProcess.cc.

186{
187 G4String thePhysicsTableFileExt;
188 if (ascii) thePhysicsTableFileExt = ".asc";
189 else thePhysicsTableFileExt = ".dat";
190
191 thePhysicsTableFileName = directory + "/";
192 thePhysicsTableFileName += tableName + "." + theProcessName + ".";
194 + thePhysicsTableFileExt;
195
197}
const G4String & GetParticleName() const
G4String thePhysicsTableFileName
Definition: G4VProcess.hh:344

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

Referenced by export_G4VProcess(), G4GammaGeneralProcess::RetrievePhysicsTable(), G4VEmProcess::RetrievePhysicsTable(), G4VEnergyLossProcess::RetrieveTable(), G4GammaGeneralProcess::StorePhysicsTable(), G4VEmProcess::StorePhysicsTable(), G4VMultipleScattering::StorePhysicsTable(), and G4VEnergyLossProcess::StoreTable().

◆ GetPILfactor()

G4double G4VProcess::GetPILfactor ( ) const
inlineinherited

Definition at line 455 of file G4VProcess.hh.

456{
457 return thePILfactor;
458}

References G4VProcess::thePILfactor.

Referenced by export_G4VProcess().

◆ GetProcessID()

size_t G4VITProcess::GetProcessID ( ) const
inlineinherited

◆ GetProcessManager()

const G4ProcessManager * G4VProcess::GetProcessManager ( )
inlinevirtualinherited

Reimplemented in G4BiasingProcessInterface, and G4WrapperProcess.

Definition at line 494 of file G4VProcess.hh.

495{
496 return aProcessManager;
497}
const G4ProcessManager * aProcessManager
Definition: G4VProcess.hh:319

References G4VProcess::aProcessManager.

Referenced by G4BiasingProcessInterface::GetProcessManager(), and G4WrapperProcess::GetProcessManager().

◆ GetProcessName()

const G4String & G4VProcess::GetProcessName ( ) const
inlineinherited

Definition at line 382 of file G4VProcess.hh.

383{
384 return theProcessName;
385}

References G4VProcess::theProcessName.

Referenced by G4VEnergyLossProcess::ActivateForcedInteraction(), G4VEmProcess::ActivateForcedInteraction(), G4ProcessManager::ActivateProcess(), G4VEmProcess::ActivateSecondaryBiasing(), G4VEnergyLossProcess::ActivateSecondaryBiasing(), G4ParallelGeometriesLimiterProcess::AddParallelWorld(), G4IonQMDPhysics::AddProcess(), G4IonINCLXXPhysics::AddProcess(), G4ProcessManager::AddProcess(), G4ProcessPlacer::AddProcessAs(), G4ITSteppingVerbose::AlongStepDoItAllDone(), G4SteppingVerbose::AlongStepDoItAllDone(), G4SteppingVerboseWithUnits::AlongStepDoItAllDone(), G4ITSteppingVerbose::AlongStepDoItOneByOne(), G4SteppingVerbose::AlongStepDoItOneByOne(), G4SteppingVerboseWithUnits::AlongStepDoItOneByOne(), G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4ITSteppingVerbose::AtRestDoItInvoked(), G4SteppingVerbose::AtRestDoItInvoked(), G4SteppingVerboseWithUnits::AtRestDoItInvoked(), G4ITSteppingVerbose::AtRestDoItOneByOne(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4HadronicProcess::BiasCrossSectionByFactor(), G4VXTRenergyLoss::BuildAngleForEnergyBank(), G4VEnergyLossProcess::BuildDEDXTable(), G4VUserPhysicsList::BuildIntegralPhysicsTable(), G4VEmProcess::BuildLambdaTable(), G4VEnergyLossProcess::BuildLambdaTable(), G4DNABrownianTransportation::BuildPhysicsTable(), G4GammaGeneralProcess::BuildPhysicsTable(), G4VEmProcess::BuildPhysicsTable(), G4VEnergyLossProcess::BuildPhysicsTable(), G4VMultipleScattering::BuildPhysicsTable(), G4LossTableManager::BuildPhysicsTable(), G4LossTableManager::BuildTables(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4ProcessManager::CheckOrderingParameters(), G4HadronicProcess::CheckResult(), G4StackChecker::ClassifyNewTrack(), G4BOptrForceCollision::ConfigureForWorker(), G4RunManagerKernel::ConfirmCoupledTransportation(), G4FastSimulationPhysics::ConstructProcess(), G4GenericBiasingPhysics::ConstructProcess(), G4IonElasticPhysics::ConstructProcess(), G4LossTableManager::CopyTables(), G4RichTrajectory::CreateAttValues(), G4RichTrajectoryPoint::CreateAttValues(), G4VPhononProcess::CreateSecondary(), G4EmExtraParameters::DefineRegParamForEM(), G4EmExtraParameters::DefineRegParamForLoss(), G4HadronicProcessStore::DeRegisterExtraProcess(), G4ITSteppingVerbose::DPSLAlongStep(), G4SteppingVerbose::DPSLAlongStep(), G4SteppingVerboseWithUnits::DPSLAlongStep(), G4ITSteppingVerbose::DPSLPostStep(), G4SteppingVerbose::DPSLPostStep(), G4SteppingVerboseWithUnits::DPSLPostStep(), G4HadronicProcessStore::Dump(), G4HadronicProcess::DumpState(), G4MuonicAtomDecay::DumpState(), G4ExceptionHandler::DumpTrackInfo(), export_G4VProcess(), G4EmCalculator::FindEmModel(), G4VEmProcess::FindLambdaMax(), G4BiasingProcessInterface::G4BiasingProcessInterface(), G4Cerenkov::G4Cerenkov(), G4ErrorEnergyLoss::G4ErrorEnergyLoss(), G4ErrorTrackLengthTarget::G4ErrorTrackLengthTarget(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4ImportanceProcess::G4ImportanceProcess(), G4MaxTimeCuts::G4MaxTimeCuts(), G4MicroElecSurface::G4MicroElecSurface(), G4MinEkineCuts::G4MinEkineCuts(), G4OpAbsorption::G4OpAbsorption(), G4OpBoundaryProcess::G4OpBoundaryProcess(), G4OpMieHG::G4OpMieHG(), G4OpRayleigh::G4OpRayleigh(), G4OpWLS::G4OpWLS(), G4OpWLS2::G4OpWLS2(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess(), G4Scintillation::G4Scintillation(), G4ScoreSplittingProcess::G4ScoreSplittingProcess(), G4SpecialCuts::G4SpecialCuts(), G4StepLimiter::G4StepLimiter(), G4UCNAbsorption::G4UCNAbsorption(), G4UCNBoundaryProcess::G4UCNBoundaryProcess(), G4UCNLoss::G4UCNLoss(), G4UCNMultiScattering::G4UCNMultiScattering(), G4UserSpecialCuts::G4UserSpecialCuts(), G4WeightCutOffProcess::G4WeightCutOffProcess(), G4WeightWindowProcess::G4WeightWindowProcess(), G4HadronicProcess::GetElementCrossSection(), G4VEmProcess::GetEmProcess(), G4GammaGeneralProcess::GetEmProcess(), G4WeightWindowProcess::GetName(), G4ProcessManager::GetProcess(), G4ProcessManager::GetProcessVectorIndex(), G4GammaGeneralProcess::GetSubProcessName(), G4ProcessManager::InActivateProcess(), G4hhIonisation::InitialiseEnergyLossProcess(), G4ProcessTable::Insert(), G4ITStepProcessor::InvokeAlongStepDoItProcs(), G4SteppingManager::InvokeAlongStepDoItProcs(), G4SteppingManager::InvokeAtRestDoItProcs(), G4SteppingManager::InvokePSDIP(), G4LossTableManager::LocalPhysicsTables(), G4ErrorPropagator::MakeOneStep(), G4VEmProcess::PostStepDoIt(), G4ITSteppingVerbose::PostStepDoItAllDone(), G4SteppingVerbose::PostStepDoItAllDone(), G4SteppingVerboseWithUnits::PostStepDoItAllDone(), G4ITSteppingVerbose::PostStepDoItOneByOne(), G4SteppingVerbose::PostStepDoItOneByOne(), G4SteppingVerboseWithUnits::PostStepDoItOneByOne(), G4VITDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4ITSteppingVerbose::PostStepVerbose(), G4EmConfigurator::PrepareModels(), G4HadronStoppingProcess::PreparePhysicsTable(), G4GammaGeneralProcess::PreparePhysicsTable(), G4VEmProcess::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), G4VMultipleScattering::PreparePhysicsTable(), G4LossTableManager::PreparePhysicsTable(), G4HadronicProcessStore::Print(), G4HadronicProcessStore::PrintHtml(), G4AnnihiToMuPair::PrintInfoDefinition(), G4GammaConversionToMuons::PrintInfoDefinition(), G4hImpactIonisation::PrintInfoDefinition(), G4ProcessPlacer::PrintProcVec(), G4VEnergyLossProcess::PrintWarning(), G4VEmProcess::PrintWarning(), G4SynchrotronRadiation::ProcessDescription(), G4Decay::ProcessDescription(), G4DecayWithSpin::ProcessDescription(), G4PionDecayMakeSpin::ProcessDescription(), G4UnknownDecay::ProcessDescription(), G4ChannelingOptrChangeCrossSection::ProposeOccurenceBiasingOperation(), G4StackManager::PushOneTrack(), G4HadronicProcessStore::Register(), G4LossTableManager::Register(), G4LossTableManager::RegisterExtraParticle(), G4HadronicProcessStore::RegisterExtraProcess(), G4HadronicProcessStore::RegisterParticle(), G4WrapperProcess::RegisterProcess(), G4PhysicsListHelper::RegisterProcess(), G4ProcessTable::Remove(), G4ParallelGeometriesLimiterProcess::RemoveParallelWorld(), G4ProcessManager::RemoveProcess(), G4ProcessPlacer::RemoveProcess(), G4GammaGeneralProcess::RetrievePhysicsTable(), G4VEmProcess::RetrievePhysicsTable(), G4VEnergyLossProcess::RetrievePhysicsTable(), G4VEmProcess::SetCrossSectionBiasingFactor(), G4VEnergyLossProcess::SetCrossSectionBiasingFactor(), G4VEnergyLossProcess::SetCSDARangeTable(), G4VEnergyLossProcess::SetInverseRangeTable(), G4VEnergyLossProcess::SetLambdaTable(), G4ProcessTableMessenger::SetNewValue(), G4ProcessTable::SetProcessActivation(), G4ProcessManager::SetProcessOrdering(), G4ProcessManager::SetProcessOrderingToFirst(), G4ProcessManager::SetProcessOrderingToLast(), G4ProcessManager::SetProcessOrderingToSecond(), G4VEnergyLossProcess::SetRangeTableForLoss(), G4VEnergyLossProcess::SetSecondaryRangeTable(), G4FastSimulationManagerProcess::SetWorldVolume(), G4ITSteppingVerbose::ShowStep(), G4SteppingVerbose::ShowStep(), G4SteppingVerboseWithUnits::ShowStep(), G4ChannelingOptrChangeCrossSection::StartRun(), G4ITSteppingVerbose::StepInfo(), G4SteppingVerbose::StepInfo(), G4SteppingVerboseWithUnits::StepInfo(), G4ITSteppingVerbose::StepInfoForLeadingTrack(), G4VEmProcess::StorePhysicsTable(), G4VMultipleScattering::StorePhysicsTable(), G4VEnergyLossProcess::StreamInfo(), G4VEmProcess::StreamInfo(), G4VMultipleScattering::StreamInfo(), G4EmCalculator::UpdateParticle(), G4ParallelWorldScoringProcess::Verbose(), G4ScoreSplittingProcess::Verbose(), G4ITSteppingVerbose::VerboseTrack(), G4SteppingVerbose::VerboseTrack(), and G4SteppingVerboseWithUnits::VerboseTrack().

◆ GetProcessState()

G4shared_ptr< G4ProcessState_Lock > G4VITProcess::GetProcessState ( )
inlineinherited

Definition at line 120 of file G4VITProcess.hh.

121 {
123 }
#define UpcastProcessState(destinationType)
Definition: G4VITProcess.hh:82

References UpcastProcessState.

◆ GetProcessSubType()

G4int G4VProcess::GetProcessSubType ( ) const
inlineinherited

◆ GetProcessType()

G4ProcessType G4VProcess::GetProcessType ( ) const
inlineinherited

◆ GetProcessTypeName()

const G4String & G4VProcess::GetProcessTypeName ( G4ProcessType  aType)
staticinherited

Definition at line 134 of file G4VProcess.cc.

135{
136 switch (aType)
137 {
138 case fNotDefined: return typeNotDefined; break;
139 case fTransportation: return typeTransportation; break;
140 case fElectromagnetic: return typeElectromagnetic; break;
141 case fOptical: return typeOptical; break;
142 case fHadronic: return typeHadronic; break;
144 case fDecay: return typeDecay; break;
145 case fGeneral: return typeGeneral; break;
146 case fParameterisation: return typeParameterisation; break;
147 case fUserDefined: return typeUserDefined; break;
148 case fPhonon: return typePhonon; break;
149 default: ;
150 }
151 return noType;
152}
@ fOptical
@ fPhonon
@ fParameterisation
@ fGeneral
@ fDecay
@ fElectromagnetic
@ fHadronic
@ fUserDefined
@ fTransportation
@ fPhotolepton_hadron
@ fNotDefined
static const G4String typeNotDefined
Definition: G4VProcess.cc:119
static const G4String typeParameterisation
Definition: G4VProcess.cc:127
static const G4String typePhotolepton_hadron
Definition: G4VProcess.cc:124
static const G4String typeElectromagnetic
Definition: G4VProcess.cc:121
static const G4String noType
Definition: G4VProcess.cc:130
static const G4String typeUserDefined
Definition: G4VProcess.cc:128
static const G4String typeDecay
Definition: G4VProcess.cc:125
static const G4String typeTransportation
Definition: G4VProcess.cc:120
static const G4String typeHadronic
Definition: G4VProcess.cc:123
static const G4String typeOptical
Definition: G4VProcess.cc:122
static const G4String typeGeneral
Definition: G4VProcess.cc:126
static const G4String typePhonon
Definition: G4VProcess.cc:129

References fDecay, fElectromagnetic, fGeneral, fHadronic, fNotDefined, fOptical, fParameterisation, fPhonon, fPhotolepton_hadron, fTransportation, fUserDefined, anonymous_namespace{G4VProcess.cc}::noType, anonymous_namespace{G4VProcess.cc}::typeDecay, anonymous_namespace{G4VProcess.cc}::typeElectromagnetic, anonymous_namespace{G4VProcess.cc}::typeGeneral, anonymous_namespace{G4VProcess.cc}::typeHadronic, anonymous_namespace{G4VProcess.cc}::typeNotDefined, anonymous_namespace{G4VProcess.cc}::typeOptical, anonymous_namespace{G4VProcess.cc}::typeParameterisation, anonymous_namespace{G4VProcess.cc}::typePhonon, anonymous_namespace{G4VProcess.cc}::typePhotolepton_hadron, anonymous_namespace{G4VProcess.cc}::typeTransportation, and anonymous_namespace{G4VProcess.cc}::typeUserDefined.

Referenced by G4RichTrajectory::CreateAttValues(), G4RichTrajectoryPoint::CreateAttValues(), G4ProcessManager::DumpInfo(), G4VProcess::DumpInfo(), G4ProcessTableMessenger::G4ProcessTableMessenger(), G4ProcessTableMessenger::GetProcessType(), G4ProcessTableMessenger::GetProcessTypeName(), and G4ProcessTableMessenger::SetNumberOfProcessType().

◆ GetState()

template<typename T >
T * G4VITProcess::GetState ( )
inlineprotectedinherited

Definition at line 212 of file G4VITProcess.hh.

213 {
214 return fpState->GetState<T>();
215 }

References G4VITProcess::fpState.

◆ GetTotalNumberOfInteractionLengthTraversed()

G4double G4VProcess::GetTotalNumberOfInteractionLengthTraversed ( ) const
inlineinherited

◆ GetVerboseLevel()

G4int G4VProcess::GetVerboseLevel ( ) const
inlineinherited

◆ InstantiateProcessState()

G4bool G4VITProcess::InstantiateProcessState ( )
inlineprotectedinherited

Definition at line 234 of file G4VITProcess.hh.

235 {
237 }
G4bool fInstantiateProcessState

References G4VITProcess::fInstantiateProcessState.

Referenced by G4VITProcess::StartTracking().

◆ isAlongStepDoItIsEnabled()

G4bool G4VProcess::isAlongStepDoItIsEnabled ( ) const
inlineinherited

Definition at line 506 of file G4VProcess.hh.

507{
508 return enableAlongStepDoIt;
509}

References G4VProcess::enableAlongStepDoIt.

Referenced by G4ProcessManager::CheckOrderingParameters().

◆ IsApplicable()

G4bool G4DNAMolecularDissociation::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 83 of file G4DNAMolecularDissociation.cc.

85{
86 if (aParticleType.GetParticleType() == "Molecule")
87 {
88#ifdef G4VERBOSE
89
90 if (fVerbose > 1)
91 {
92 G4cout << "G4MolecularDissociation::IsApplicable(";
93 G4cout << aParticleType.GetParticleName() << ",";
94 G4cout << aParticleType.GetParticleType() << ")" << G4endl;
95 }
96#endif
97 return (true);
98 }
99 else
100 {
101 return false;
102 }
103}
const G4String & GetParticleType() const

References fVerbose, G4cout, G4endl, G4ParticleDefinition::GetParticleName(), and G4ParticleDefinition::GetParticleType().

◆ isAtRestDoItIsEnabled()

G4bool G4VProcess::isAtRestDoItIsEnabled ( ) const
inlineinherited

Definition at line 500 of file G4VProcess.hh.

501{
502 return enableAtRestDoIt;
503}

References G4VProcess::enableAtRestDoIt.

Referenced by G4ProcessManager::CheckOrderingParameters().

◆ isPostStepDoItIsEnabled()

G4bool G4VProcess::isPostStepDoItIsEnabled ( ) const
inlineinherited

Definition at line 512 of file G4VProcess.hh.

513{
514 return enablePostStepDoIt;
515}

References G4VProcess::enablePostStepDoIt.

Referenced by G4ProcessManager::CheckOrderingParameters().

◆ operator!=() [1/2]

G4bool G4VITProcess::operator!= ( const G4VITProcess right) const
inherited

◆ operator!=() [2/2]

G4bool G4VProcess::operator!= ( const G4VProcess right) const
inherited

Definition at line 161 of file G4VProcess.cc.

162{
163 return (this != &right);
164}

◆ operator=()

G4DNAMolecularDissociation & G4DNAMolecularDissociation::operator= ( const G4DNAMolecularDissociation right)
delete

◆ operator==() [1/2]

G4bool G4VITProcess::operator== ( const G4VITProcess right) const
inherited

◆ operator==() [2/2]

G4bool G4VProcess::operator== ( const G4VProcess right) const
inherited

Definition at line 155 of file G4VProcess.cc.

156{
157 return (this == &right);
158}

◆ PostStepDoIt()

G4VParticleChange * G4DNAMolecularDissociation::PostStepDoIt ( const G4Track track,
const G4Step step 
)
overridevirtual

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 338 of file G4DNAMolecularDissociation.cc.

340{
341 return AtRestDoIt(track, step);
342}
G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &step) override

References AtRestDoIt().

◆ PostStepGetPhysicalInteractionLength()

G4double G4DNAMolecularDissociation::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
overridevirtual

Reimplemented from G4VITRestDiscreteProcess.

Definition at line 299 of file G4DNAMolecularDissociation.cc.

302{
303 return 0; //-1*(-log(1-G4UniformRand())*100*1e-15*s);
304}

◆ PostStepGPIL()

G4double G4VProcess::PostStepGPIL ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
inlineinherited

Definition at line 479 of file G4VProcess.hh.

482{
483 return thePILfactor *
484 PostStepGetPhysicalInteractionLength(track, previousStepSize, condition);
485}
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0

References condition(), G4VProcess::PostStepGetPhysicalInteractionLength(), and G4VProcess::thePILfactor.

Referenced by G4SteppingManager::DefinePhysicalStepLength(), and G4ITStepProcessor::DoDefinePhysicalStepLength().

◆ PreparePhysicsTable()

virtual void G4VProcess::PreparePhysicsTable ( const G4ParticleDefinition )
inlinevirtualinherited

◆ PrepareWorkerPhysicsTable()

void G4VProcess::PrepareWorkerPhysicsTable ( const G4ParticleDefinition part)
virtualinherited

Reimplemented in G4BiasingProcessInterface.

Definition at line 206 of file G4VProcess.cc.

207{
209}
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:194

References G4VProcess::PreparePhysicsTable().

Referenced by G4BiasingProcessInterface::PrepareWorkerPhysicsTable().

◆ ProcessDescription()

void G4VProcess::ProcessDescription ( std::ostream &  outfile) const
virtualinherited

Reimplemented in G4AdjointAlongStepWeightCorrection, G4AdjointForcedInteractionForGamma, G4AdjointhMultipleScattering, G4ContinuousGainOfEnergy, G4eAdjointMultipleScattering, G4eInverseBremsstrahlung, G4eInverseCompton, G4eInverseIonisation, G4InversePEEffect, G4IonInverseIonisation, G4eeToHadrons, G4hBremsstrahlung, G4hhIonisation, G4hPairProduction, G4mplIonisation, G4RayleighScattering, G4ePairProduction, G4MuBremsstrahlung, G4MuIonisation, G4MuMultipleScattering, G4MuPairProduction, G4PolarizedAnnihilation, G4PolarizedBremsstrahlung, G4PolarizedCompton, G4PolarizedGammaConversion, G4PolarizedIonisation, G4PolarizedPhotoElectric, G4ComptonScattering, G4CoulombScattering, G4eBremsstrahlung, G4eIonisation, G4eMultipleScattering, G4eplusAnnihilation, G4GammaConversion, G4hIonisation, G4hMultipleScattering, G4ionIonisation, G4NuclearStopping, G4PhotoElectricEffect, G4ForwardXrayTR, G4GammaXTRadiator, G4GaussXTRadiator, G4RegularXTRadiator, G4Scintillation, G4StrawTubeXTRadiator, G4SynchrotronRadiation, G4TransitionRadiation, G4TransparentRegXTRadiator, G4VTransitionRadiation, G4VXTRenergyLoss, G4XTRGammaRadModel, G4XTRRegularRadModel, G4XTRTransparentRegRadModel, G4Cerenkov, G4Radioactivation, G4RadioactiveDecay, G4ElectronNuclearProcess, G4MuonNuclearProcess, G4NeutrinoElectronProcess, G4NeutronFissionProcess, G4PositronNuclearProcess, G4HadronicAbsorptionBertini, G4HadronicAbsorptionFritiof, G4HadronicAbsorptionFritiofWithBinaryCascade, G4HadronStoppingProcess, G4MuonicAtomDecay, G4MuonMinusAtomicCapture, G4MuonMinusCapture, G4Transportation, G4NeutronCaptureProcess, G4GammaGeneralProcess, G4Decay, G4DecayWithSpin, G4PionDecayMakeSpin, G4UnknownDecay, G4VEmProcess, G4VEnergyLossProcess, G4VMultipleScattering, G4HadronicProcess, G4ElNeutrinoNucleusProcess, G4HadronElasticProcess, and G4MuNeutrinoNucleusProcess.

Definition at line 175 of file G4VProcess.cc.

176{
177 outFile << "This process has not yet been described\n";
178}

Referenced by G4LossTableManager::DumpHtml(), G4HadronicProcessStore::PrintHtml(), and G4GammaGeneralProcess::ProcessDescription().

◆ ProposesTimeStep()

G4bool G4VITProcess::ProposesTimeStep ( ) const
inlineinherited

◆ ResetNumberOfInteractionLengthLeft()

void G4VITProcess::ResetNumberOfInteractionLengthLeft ( )
inlinevirtualinherited

◆ ResetProcessState()

void G4VITProcess::ResetProcessState ( )
inlineinherited

◆ RetrievePhysicsTable()

virtual G4bool G4VProcess::RetrievePhysicsTable ( const G4ParticleDefinition ,
const G4String ,
G4bool   
)
inlinevirtualinherited

◆ RetrieveProcessInfo()

void G4VITProcess::RetrieveProcessInfo ( )
protectedinherited

◆ SetDisplacer()

void G4DNAMolecularDissociation::SetDisplacer ( Species pSpecies,
Displacer pDisplacer 
)

◆ SetInstantiateProcessState()

void G4VITProcess::SetInstantiateProcessState ( G4bool  flag)
inlineprotectedinherited

◆ SetMasterProcess()

void G4VProcess::SetMasterProcess ( G4VProcess masterP)
virtualinherited

◆ SetPILfactor()

void G4VProcess::SetPILfactor ( G4double  value)
inlineinherited

Definition at line 449 of file G4VProcess.hh.

450{
451 if (value>0.) { thePILfactor = value; }
452}

References G4VProcess::thePILfactor.

Referenced by export_G4VProcess().

◆ SetProcessManager()

void G4VProcess::SetProcessManager ( const G4ProcessManager procMan)
inlinevirtualinherited

◆ SetProcessState()

void G4VITProcess::SetProcessState ( G4shared_ptr< G4ProcessState_Lock aProcInfo)
inlineinherited

◆ SetProcessSubType()

void G4VProcess::SetProcessSubType ( G4int  value)
inlineinherited

Definition at line 406 of file G4VProcess.hh.

407{
408 theProcessSubType = value;
409}

References G4VProcess::theProcessSubType.

Referenced by G4DNAElectronHoleRecombination::Create(), G4DNASecondOrderReaction::Create(), G4AnnihiToMuPair::G4AnnihiToMuPair(), G4BiasingProcessInterface::G4BiasingProcessInterface(), G4Cerenkov::G4Cerenkov(), G4ComptonScattering::G4ComptonScattering(), G4CoulombScattering::G4CoulombScattering(), G4CoupledTransportation::G4CoupledTransportation(), G4Decay::G4Decay(), G4DecayWithSpin::G4DecayWithSpin(), G4DNAAttachment::G4DNAAttachment(), G4DNABrownianTransportation::G4DNABrownianTransportation(), G4DNAChargeDecrease::G4DNAChargeDecrease(), G4DNAChargeIncrease::G4DNAChargeIncrease(), G4DNAElastic::G4DNAElastic(), G4DNAElectronSolvation::G4DNAElectronSolvation(), G4DNAExcitation::G4DNAExcitation(), G4DNAIonisation::G4DNAIonisation(), G4DNAMolecularDissociation(), G4DNAScavengerProcess::G4DNAScavengerProcess(), G4DNAVibExcitation::G4DNAVibExcitation(), G4eBremsstrahlung::G4eBremsstrahlung(), G4eeToHadrons::G4eeToHadrons(), G4eIonisation::G4eIonisation(), G4ePairProduction::G4ePairProduction(), G4eplusAnnihilation::G4eplusAnnihilation(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4GammaConversion::G4GammaConversion(), G4GammaConversionToMuons::G4GammaConversionToMuons(), G4GammaGeneralProcess::G4GammaGeneralProcess(), G4HadronicProcess::G4HadronicProcess(), G4hhIonisation::G4hhIonisation(), G4hIonisation::G4hIonisation(), G4ionIonisation::G4ionIonisation(), G4ITTransportation::G4ITTransportation(), G4JAEAElasticScattering::G4JAEAElasticScattering(), G4MicroElecElastic::G4MicroElecElastic(), G4MicroElecInelastic::G4MicroElecInelastic(), G4MicroElecLOPhononScattering::G4MicroElecLOPhononScattering(), G4MicroElecSurface::G4MicroElecSurface(), G4mplIonisation::G4mplIonisation(), G4MuBremsstrahlung::G4MuBremsstrahlung(), G4MuIonisation::G4MuIonisation(), G4MuonMinusAtomicCapture::G4MuonMinusAtomicCapture(), G4MuPairProduction::G4MuPairProduction(), G4NeutronKiller::G4NeutronKiller(), G4NuclearStopping::G4NuclearStopping(), G4OpAbsorption::G4OpAbsorption(), G4OpBoundaryProcess::G4OpBoundaryProcess(), G4OpMieHG::G4OpMieHG(), G4OpRayleigh::G4OpRayleigh(), G4OpWLS::G4OpWLS(), G4OpWLS2::G4OpWLS2(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4PhotoElectricEffect::G4PhotoElectricEffect(), G4PionDecayMakeSpin::G4PionDecayMakeSpin(), G4PolarizedCompton::G4PolarizedCompton(), G4PolarizedGammaConversion::G4PolarizedGammaConversion(), G4PolarizedIonisation::G4PolarizedIonisation(), G4PolarizedPhotoElectric::G4PolarizedPhotoElectric(), G4RadioactiveDecay::G4RadioactiveDecay(), G4RayleighScattering::G4RayleighScattering(), G4Scintillation::G4Scintillation(), G4StepLimiter::G4StepLimiter(), G4SynchrotronRadiation::G4SynchrotronRadiation(), G4SynchrotronRadiationInMat::G4SynchrotronRadiationInMat(), G4TransitionRadiation::G4TransitionRadiation(), G4Transportation::G4Transportation(), G4UCNAbsorption::G4UCNAbsorption(), G4UCNBoundaryProcess::G4UCNBoundaryProcess(), G4UCNLoss::G4UCNLoss(), G4UCNMultiScattering::G4UCNMultiScattering(), G4UnknownDecay::G4UnknownDecay(), G4UserSpecialCuts::G4UserSpecialCuts(), G4VMultipleScattering::G4VMultipleScattering(), G4VTransitionRadiation::G4VTransitionRadiation(), G4VXTRenergyLoss::G4VXTRenergyLoss(), and G4Decay::SetExtDecayer().

◆ SetProcessType()

void G4VProcess::SetProcessType ( G4ProcessType  aType)
inlineinherited

Definition at line 394 of file G4VProcess.hh.

395{
396 theProcessType = aType;
397}

References G4VProcess::theProcessType.

Referenced by G4MaxTimeCuts::G4MaxTimeCuts(), and G4MinEkineCuts::G4MinEkineCuts().

◆ SetVerbose()

void G4DNAMolecularDissociation::SetVerbose ( G4int  verbose)

Definition at line 308 of file G4DNAMolecularDissociation.cc.

309{
310 fVerbose = verbose;
311}

References fVerbose.

◆ SetVerboseLevel()

void G4VProcess::SetVerboseLevel ( G4int  value)
inlineinherited

◆ StartTracking()

void G4VITProcess::StartTracking ( G4Track track)
virtualinherited

Reimplemented from G4VProcess.

Reimplemented in G4DNASecondOrderReaction, G4DNAElectronHoleRecombination, G4DNAScavengerProcess, G4ITTransportation, and G4DNABrownianTransportation.

Definition at line 85 of file G4VITProcess.cc.

86{
87 G4TrackingInformation* trackingInfo = GetIT(track)->GetTrackingInfo();
89 {
90 // fpState = new G4ProcessState();
91 fpState.reset(new G4ProcessState());
92 }
93
94 theNumberOfInteractionLengthLeft = &(fpState->theNumberOfInteractionLengthLeft );
95 theInteractionTimeLeft = &(fpState->theInteractionTimeLeft );
96 currentInteractionLength = &(fpState->currentInteractionLength );
98 // fpState = 0;
99 fpState.reset();
100}
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
G4TrackingInformation * GetTrackingInfo()
Definition: G4IT.hh:143
void RecordProcessState(G4shared_ptr< G4ProcessState_Lock >, size_t index)
G4double * theInteractionTimeLeft
G4double * theNumberOfInteractionLengthLeft
G4bool InstantiateProcessState()
G4double * currentInteractionLength

References G4VITProcess::currentInteractionLength, G4VITProcess::fProcessID, G4VITProcess::fpState, GetIT(), G4IT::GetTrackingInfo(), G4VITProcess::InstantiateProcessState(), G4TrackingInformation::RecordProcessState(), G4VITProcess::theInteractionTimeLeft, and G4VITProcess::theNumberOfInteractionLengthLeft.

Referenced by G4DNASecondOrderReaction::StartTracking(), G4DNAElectronHoleRecombination::StartTracking(), G4DNAScavengerProcess::StartTracking(), and G4ITTransportation::StartTracking().

◆ StorePhysicsTable()

virtual G4bool G4VProcess::StorePhysicsTable ( const G4ParticleDefinition ,
const G4String ,
G4bool   
)
inlinevirtualinherited

◆ SubtractNumberOfInteractionLengthLeft()

void G4VITProcess::SubtractNumberOfInteractionLengthLeft ( G4double  previousStepSize)
inlineprotectedvirtualinherited

Definition at line 292 of file G4VITProcess.hh.

293{
294 if (fpState->currentInteractionLength > 0.0)
295 {
296 fpState->theNumberOfInteractionLengthLeft -= previousStepSize
297 / fpState->currentInteractionLength;
298 if (fpState->theNumberOfInteractionLengthLeft < 0.)
299 {
300 fpState->theNumberOfInteractionLengthLeft = CLHEP::perMillion;
301 }
302
303 }
304 else
305 {
306#ifdef G4VERBOSE
307 if (verboseLevel > 0)
308 {
309 G4cerr << "G4VITProcess::SubtractNumberOfInteractionLengthLeft()";
310 G4cerr << " [" << theProcessName << "]" << G4endl;
311 G4cerr << " currentInteractionLength = "
312 << fpState->currentInteractionLength << " [mm]";
313 G4cerr << " previousStepSize = " << previousStepSize << " [mm]";
314 G4cerr << G4endl;
315 }
316#endif
317 G4String msg = "Negative currentInteractionLength for ";
318 msg += theProcessName;
319 G4Exception("G4VITProcess::SubtractNumberOfInteractionLengthLeft()",
320 "ProcMan201",EventMustBeAborted,
321 msg);
322 }
323}
@ EventMustBeAborted
G4GLOB_DLL std::ostream G4cerr
static constexpr double perMillion

References EventMustBeAborted, G4VITProcess::fpState, G4cerr, G4endl, G4Exception(), CLHEP::perMillion, G4VProcess::theProcessName, and G4VProcess::verboseLevel.

Referenced by G4VITDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), and G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength().

Field Documentation

◆ aParticleChange

G4ParticleChange G4VProcess::aParticleChange
protectedinherited

◆ aProcessManager

const G4ProcessManager* G4VProcess::aProcessManager = nullptr
protectedinherited

Definition at line 319 of file G4VProcess.hh.

Referenced by G4VProcess::GetProcessManager(), and G4VProcess::SetProcessManager().

◆ currentInteractionLength

G4double* G4VITProcess::currentInteractionLength
privateinherited

Definition at line 254 of file G4VITProcess.hh.

Referenced by G4VITProcess::G4VITProcess(), and G4VITProcess::StartTracking().

◆ enableAlongStepDoIt

G4bool G4VProcess::enableAlongStepDoIt = true
protectedinherited

◆ enableAtRestDoIt

G4bool G4VProcess::enableAtRestDoIt = true
protectedinherited

◆ enablePostStepDoIt

G4bool G4VProcess::enablePostStepDoIt = true
protectedinherited

◆ fDecayAtFixedTime

G4bool G4DNAMolecularDissociation::fDecayAtFixedTime
private

◆ fDisplacementMap

DisplacementMap G4DNAMolecularDissociation::fDisplacementMap
private

Definition at line 107 of file G4DNAMolecularDissociation.hh.

Referenced by DecayIt(), GetDisplacer(), and SetDisplacer().

◆ fInstantiateProcessState

G4bool G4VITProcess::fInstantiateProcessState
privateinherited

◆ fNbProcess

size_t * G4VITProcess::fNbProcess = 0
privateinherited

◆ fProcessID

size_t G4VITProcess::fProcessID
privateinherited

◆ fProcessTable

G4ProcessTable* G4VProcess::fProcessTable = nullptr
privateinherited

Definition at line 374 of file G4VProcess.hh.

Referenced by G4VProcess::G4VProcess(), and G4VProcess::~G4VProcess().

◆ fProposesTimeStep

G4bool G4VITProcess::fProposesTimeStep
protectedinherited

◆ fpState

G4shared_ptr<G4ProcessState> G4VITProcess::fpState
protectedinherited

◆ fVerbose

G4int G4DNAMolecularDissociation::fVerbose
private

◆ masterProcessShadow

G4VProcess* G4VProcess::masterProcessShadow = nullptr
privateinherited

Definition at line 370 of file G4VProcess.hh.

Referenced by G4VProcess::GetMasterProcess(), and G4VProcess::SetMasterProcess().

◆ pParticleChange

G4VParticleChange* G4VProcess::pParticleChange = nullptr
protectedinherited

Definition at line 321 of file G4VProcess.hh.

Referenced by G4VMultipleScattering::AddEmModel(), G4VEmProcess::AddEmModel(), G4VEnergyLossProcess::AddEmModel(), G4ImportanceProcess::AlongStepDoIt(), G4WeightCutOffProcess::AlongStepDoIt(), G4WeightWindowProcess::AlongStepDoIt(), G4VContinuousDiscreteProcess::AlongStepDoIt(), G4VContinuousProcess::AlongStepDoIt(), G4VRestContinuousDiscreteProcess::AlongStepDoIt(), G4VRestContinuousProcess::AlongStepDoIt(), G4ParallelWorldProcess::AlongStepDoIt(), G4ParallelWorldScoringProcess::AlongStepDoIt(), G4VITRestProcess::AtRestDoIt(), G4VRestContinuousDiscreteProcess::AtRestDoIt(), G4VRestContinuousProcess::AtRestDoIt(), G4VRestDiscreteProcess::AtRestDoIt(), G4VRestProcess::AtRestDoIt(), G4ParallelWorldProcess::AtRestDoIt(), G4ParallelWorldScoringProcess::AtRestDoIt(), G4ScoreSplittingProcess::AtRestDoIt(), G4VITRestDiscreteProcess::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4DNAElectronHoleRecombination::Create(), G4DNASecondOrderReaction::Create(), G4VEnergyLossProcess::FillSecondariesAlongStep(), G4Decay::G4Decay(), G4DNAMolecularDissociation(), G4DNAScavengerProcess::G4DNAScavengerProcess(), G4ImportanceProcess::G4ImportanceProcess(), G4ITTransportation::G4ITTransportation(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess(), G4RadioactiveDecay::G4RadioactiveDecay(), G4ScoreSplittingProcess::G4ScoreSplittingProcess(), G4Transportation::G4Transportation(), G4UnknownDecay::G4UnknownDecay(), G4VEmProcess::G4VEmProcess(), G4VEnergyLossProcess::G4VEnergyLossProcess(), G4VMultipleScattering::G4VMultipleScattering(), G4VProcess::G4VProcess(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4WeightCutOffProcess::G4WeightCutOffProcess(), G4WeightWindowProcess::G4WeightWindowProcess(), G4VITDiscreteProcess::PostStepDoIt(), G4VContinuousDiscreteProcess::PostStepDoIt(), G4VDiscreteProcess::PostStepDoIt(), G4VRestContinuousDiscreteProcess::PostStepDoIt(), G4VRestDiscreteProcess::PostStepDoIt(), G4ParallelWorldProcess::PostStepDoIt(), G4ParallelWorldScoringProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4NeutronKiller::PostStepDoIt(), G4VITRestDiscreteProcess::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), and G4VTransitionRadiation::PostStepDoIt().

◆ theInitialNumberOfInteractionLength

G4double G4VProcess::theInitialNumberOfInteractionLength = -1.0
protectedinherited

◆ theInteractionTimeLeft

G4double* G4VITProcess::theInteractionTimeLeft
privateinherited

◆ theNumberOfInteractionLengthLeft

G4double* G4VITProcess::theNumberOfInteractionLengthLeft
privateinherited

Definition at line 253 of file G4VITProcess.hh.

Referenced by G4VITProcess::G4VITProcess(), and G4VITProcess::StartTracking().

◆ thePhysicsTableFileName

G4String G4VProcess::thePhysicsTableFileName
protectedinherited

Definition at line 344 of file G4VProcess.hh.

Referenced by G4VProcess::GetPhysicsTableFileName().

◆ thePILfactor

G4double G4VProcess::thePILfactor = 1.0
protectedinherited

◆ theProcessName

G4String G4VProcess::theProcessName
protectedinherited

◆ theProcessSubType

G4int G4VProcess::theProcessSubType = -1
protectedinherited

◆ theProcessType

G4ProcessType G4VProcess::theProcessType = fNotDefined
protectedinherited

◆ verboseLevel

G4int G4VProcess::verboseLevel = 0
protectedinherited

Definition at line 356 of file G4VProcess.hh.

Referenced by G4VEnergyLossProcess::ActivateForcedInteraction(), G4VEmProcess::ActivateForcedInteraction(), G4VEmProcess::ActivateSecondaryBiasing(), G4VEnergyLossProcess::ActivateSecondaryBiasing(), G4LowECapture::AddRegion(), G4CoupledTransportation::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4hImpactIonisation::AntiProtonParametrisedDEDX(), G4ParallelWorldScoringProcess::AtRestDoIt(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4VXTRenergyLoss::BuildAngleForEnergyBank(), G4VXTRenergyLoss::BuildAngleTable(), G4VEnergyLossProcess::BuildDEDXTable(), G4VXTRenergyLoss::BuildEnergyTable(), G4VXTRenergyLoss::BuildGlobalAngleTable(), G4VEmProcess::BuildLambdaTable(), G4hImpactIonisation::BuildLambdaTable(), G4VEnergyLossProcess::BuildLambdaTable(), G4hImpactIonisation::BuildLossTable(), G4DNABrownianTransportation::BuildPhysicsTable(), G4GammaGeneralProcess::BuildPhysicsTable(), G4LowECapture::BuildPhysicsTable(), G4VEmProcess::BuildPhysicsTable(), G4VEnergyLossProcess::BuildPhysicsTable(), G4VMultipleScattering::BuildPhysicsTable(), G4SynchrotronRadiation::BuildPhysicsTable(), G4VXTRenergyLoss::BuildPhysicsTable(), G4hImpactIonisation::BuildPhysicsTable(), G4ChargeExchangeProcess::BuildPhysicsTable(), G4OpRayleigh::CalculateRayleighMeanFreePaths(), G4PolarizedAnnihilation::ComputeSaturationFactor(), G4PolarizedCompton::ComputeSaturationFactor(), G4PolarizedIonisation::ComputeSaturationFactor(), G4DNAElectronHoleRecombination::Create(), G4DNASecondOrderReaction::Create(), G4VPhononProcess::CreateSecondary(), G4VProcess::EndTracking(), G4VEmProcess::FindLambdaMax(), G4Cerenkov::G4Cerenkov(), G4ChargeExchangeProcess::G4ChargeExchangeProcess(), G4CoupledTransportation::G4CoupledTransportation(), G4DNAMolecularDissociation(), G4DNAScavengerProcess::G4DNAScavengerProcess(), G4ErrorEnergyLoss::G4ErrorEnergyLoss(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4GaussXTRadiator::G4GaussXTRadiator(), G4ImportanceProcess::G4ImportanceProcess(), G4MaxTimeCuts::G4MaxTimeCuts(), G4MicroElecSurface::G4MicroElecSurface(), G4MinEkineCuts::G4MinEkineCuts(), G4OpAbsorption::G4OpAbsorption(), G4OpBoundaryProcess::G4OpBoundaryProcess(), G4OpMieHG::G4OpMieHG(), G4OpRayleigh::G4OpRayleigh(), G4OpWLS::G4OpWLS(), G4OpWLS2::G4OpWLS2(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess(), G4PolarizedIonisation::G4PolarizedIonisation(), G4Scintillation::G4Scintillation(), G4ScoreSplittingProcess::G4ScoreSplittingProcess(), G4SpecialCuts::G4SpecialCuts(), G4StepLimiter::G4StepLimiter(), G4StrawTubeXTRadiator::G4StrawTubeXTRadiator(), G4SynchrotronRadiation::G4SynchrotronRadiation(), G4TransparentRegXTRadiator::G4TransparentRegXTRadiator(), G4Transportation::G4Transportation(), G4UCNAbsorption::G4UCNAbsorption(), G4UCNBoundaryProcess::G4UCNBoundaryProcess(), G4UCNLoss::G4UCNLoss(), G4UCNMultiScattering::G4UCNMultiScattering(), G4UserSpecialCuts::G4UserSpecialCuts(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4WeightCutOffProcess::G4WeightCutOffProcess(), G4WeightWindowProcess::G4WeightWindowProcess(), G4VXTRenergyLoss::GetAngleVector(), G4ChargeExchangeProcess::GetElementCrossSection(), G4VXTRenergyLoss::GetGasZmuProduct(), G4PhononDownconversion::GetMeanFreePath(), G4PhononScattering::GetMeanFreePath(), G4PolarizedCompton::GetMeanFreePath(), G4VXTRenergyLoss::GetMeanFreePath(), G4UCNAbsorption::GetMeanFreePath(), G4PolarizedAnnihilation::GetMeanFreePath(), G4PolarizedIonisation::GetMeanFreePath(), G4SynchrotronRadiation::GetMeanFreePath(), G4VXTRenergyLoss::GetNumberOfPhotons(), G4VXTRenergyLoss::GetPlateZmuProduct(), G4SynchrotronRadiation::GetRandomEnergySR(), G4VProcess::GetVerboseLevel(), G4hhIonisation::InitialiseEnergyLossProcess(), G4eeToHadrons::InitialiseProcess(), G4hImpactIonisation::InitializeMe(), G4UCNBoundaryProcess::MRreflect(), G4UCNBoundaryProcess::MRreflectHigh(), G4DNASecondOrderReaction::PostStepDoIt(), G4ParallelWorldScoringProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4DNAScavengerProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNAbsorption::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UCNLoss::PostStepDoIt(), G4UCNMultiScattering::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpMieHG::PostStepDoIt(), G4OpRayleigh::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4VITDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength(), G4PolarizedAnnihilation::PostStepGetPhysicalInteractionLength(), G4PolarizedCompton::PostStepGetPhysicalInteractionLength(), G4PolarizedIonisation::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4GammaGeneralProcess::PreparePhysicsTable(), G4VEmProcess::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), G4VMultipleScattering::PreparePhysicsTable(), G4hImpactIonisation::ProtonParametrisedDEDX(), G4UCNBoundaryProcess::Reflect(), G4CoupledTransportation::ReportInexactEnergy(), G4CoupledTransportation::ReportMissingLogger(), G4GammaGeneralProcess::RetrievePhysicsTable(), G4VEmProcess::RetrievePhysicsTable(), G4VEnergyLossProcess::RetrievePhysicsTable(), G4VEnergyLossProcess::RetrieveTable(), G4VEmProcess::SetCrossSectionBiasingFactor(), G4VEnergyLossProcess::SetCrossSectionBiasingFactor(), G4VEnergyLossProcess::SetCSDARangeTable(), G4CoupledTransportation::SetHighLooperThresholds(), G4Transportation::SetHighLooperThresholds(), G4VEnergyLossProcess::SetInverseRangeTable(), G4LowECapture::SetKinEnergyLimit(), G4NeutronKiller::SetKinEnergyLimit(), G4VEnergyLossProcess::SetLambdaTable(), G4CoupledTransportation::SetLowLooperThresholds(), G4Transportation::SetLowLooperThresholds(), G4VEnergyLossProcess::SetRangeTableForLoss(), G4VEnergyLossProcess::SetSecondaryRangeTable(), G4NeutronKiller::SetTimeLimit(), G4VProcess::SetVerboseLevel(), G4Cerenkov::SetVerboseLevel(), G4Scintillation::SetVerboseLevel(), G4OpAbsorption::SetVerboseLevel(), G4OpBoundaryProcess::SetVerboseLevel(), G4OpMieHG::SetVerboseLevel(), G4OpRayleigh::SetVerboseLevel(), G4OpWLS::SetVerboseLevel(), G4OpWLS2::SetVerboseLevel(), G4FastSimulationManagerProcess::SetWorldVolume(), G4GaussXTRadiator::SpectralXTRdEdx(), G4RegularXTRadiator::SpectralXTRdEdx(), G4TransparentRegXTRadiator::SpectralXTRdEdx(), G4XTRRegularRadModel::SpectralXTRdEdx(), G4VProcess::StartTracking(), G4CoupledTransportation::StartTracking(), G4VEmProcess::StorePhysicsTable(), G4VMultipleScattering::StorePhysicsTable(), G4VEnergyLossProcess::StoreTable(), G4VEnergyLossProcess::StreamInfo(), G4VEmProcess::StreamInfo(), G4VMultipleScattering::StreamInfo(), G4VITProcess::SubtractNumberOfInteractionLengthLeft(), and G4VProcess::SubtractNumberOfInteractionLengthLeft().


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