Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VITRestDiscreteProcess.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: G4VITRestDiscreteProcess.hh 64057 2012-10-30 15:04:49Z gcosmo $
27//
28/// \brief Identical to G4VRestDiscreteProcess with dependency from G4VITProcess
29//
30// WARNING : This class is released as a prototype.
31// It might strongly evolve or even disapear in the next releases.
32//
33// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
34//
35// History:
36// -----------
37// 10 Oct 2011 M.Karamitros created
38//
39// -------------------------------------------------------------------
40
41#ifndef G4VITRESTDISCRETEPROCESS_H
42#define G4VITRESTDISCRETEPROCESS_H
43
45
46#include "G4VITProcess.hh"
47
48/**
49 * Identical to G4VRestDiscreteProcess with dependency from G4VITProcess
50 */
51
53{
54 // Abstract class which defines the public behavior of
55 // rest + discrete physics interactions.
56 public:
57
59 G4ProcessType aType = fNotDefined );
61
63
64 public :// with description
66 const G4Track& track,
67 G4double previousStepSize,
69 );
70
72 const G4Track& ,
73 const G4Step&
74 );
75
77 const G4Track& ,
79 );
80
82 const G4Track& ,
83 const G4Step&
84 );
85
86 // no operation in AlongStepDoIt
88 const G4Track&,
89 G4double ,
90 G4double ,
91 G4double& ,
93 ){ return -1.0; }
94
95 // no operation in AlongStepDoIt
97 const G4Track& ,
98 const G4Step&
99 ) {return 0;}
100
101 protected:// with description
102 virtual G4double GetMeanFreePath(const G4Track& aTrack,
103 G4double previousStepSize,
105 )=0;
106 // Calculates from the macroscopic cross section a mean
107 // free path, the value is returned in units of distance.
108
110 // Calculates the mean life-time (i.e. for decays) of the
111 // particle at rest due to the occurence of the given process,
112 // or converts the probability of interaction (i.e. for
113 // annihilation) into the life-time of the particle for the
114 // occurence of the given process.
115
116 private:
117 // hide default constructor and assignment operator as private
119 G4VITRestDiscreteProcess & operator=(const G4VITRestDiscreteProcess &right);
120
121};
122
123// -----------------------------------------
124// inlined function members implementation
125// -----------------------------------------
126inline
128 const G4Track& track,
129 G4double previousStepSize,
131 )
132{
133 if ( (previousStepSize < 0.0) || (fpState->theNumberOfInteractionLengthLeft<=0.0)) {
134 // beggining of tracking (or just after DoIt of this process)
136 } else if ( previousStepSize > 0.0) {
137 // subtract NumberOfInteractionLengthLeft
139 } else {
140 // zero step
141 // DO NOTHING
142 }
143
144 // condition is set to "Not Forced"
146
147 // get mean free path
148 fpState->currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
149
150 G4double value;
153 } else {
154 value = DBL_MAX;
155 }
156#ifdef G4VERBOSE
157 if (verboseLevel>1){
158 G4cout << "G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength ";
159 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
160 track.GetDynamicParticle()->DumpInfo();
161 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
162 G4cout << "InteractionLength= " << value/CLHEP::cm <<"[cm] " <<G4endl;
163 }
164#endif
165 return value;
166}
167
169 const G4Track& ,
170 const G4Step&
171 )
172{
173// reset NumberOfInteractionLengthLeft
175
176 return pParticleChange;
177}
178
180 const G4Track& track,
182 )
183{
184 // beggining of tracking
186
187 // condition is set to "Not Forced"
189
190 // get mean life time
192
193#ifdef G4VERBOSE
195 G4cout << "G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength ";
196 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
197 track.GetDynamicParticle()->DumpInfo();
198 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
199 G4cout << "MeanLifeTime = " << fpState->currentInteractionLength/CLHEP::ns << "[ns]" <<G4endl;
200 }
201#endif
202
204}
205
206
208 const G4Track&,
209 const G4Step&
210 )
211{
212// clear NumberOfInteractionLengthLeft
214
215 return pParticleChange;
216}
217
218
219
220#endif
221
222
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
G4GPILSelection
G4ProcessType
@ fNotDefined
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void DumpInfo(G4int mode=0) const
const G4String & GetName() const
Definition: G4Material.hh:177
Definition: G4Step.hh:78
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
Definition: G4VITProcess.cc:96
G4ProcessState * fpState
virtual void ResetNumberOfInteractionLengthLeft()
virtual void ClearNumberOfInteractionLengthLeft()
Identical to G4VRestDiscreteProcess with dependency from G4VITProcess.
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)=0
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
#define DBL_MAX
Definition: templates.hh:83