Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4IVContinuousDiscreteProcess.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//
27// $Id$
28//
29// $Id:
30// ------------------------------------------------------------
31// GEANT 4 class header file
32//
33//
34// Class Description
35// Abstract class which defines the public behavior of
36// discrete physics interactions using integral approach
37//
38// ------------------------------------------------------------
39// New Physics scheme 8 Mar. 1997 H.Kurahige
40// ------------------------------------------------------------
41// modified 26 Mar. 1997 H.Kurashige
42// modified 16 Apr. 1997 L.Urban
43// modified AlongStepGPIL etc. 17 Dec. 1997 H.Kurashige
44// fix bugs in GetGPILSelection() 24 Jan. 1998 H.Kurashige
45// modified for new ParticleChange 12 Mar. 1998 H.Kurashige
46
47#ifndef G4IVContinuousDiscreteProcess_h
48#define G4IVContinuousDiscreteProcess_h 1
49
51
52#include "globals.hh"
53#include "G4ios.hh"
54
55#include "G4VProcess.hh"
56#include "G4Track.hh"
57#include "G4Step.hh"
58#include "G4MaterialTable.hh"
59
61{
62 // Abstract class which defines the public behavior of
63 // discrete physics interactions.
64 public:
65
67 G4ProcessType aType = fNotDefined );
69
71
72 public: // with description
74 const G4Track& track,
75 G4double previousStepSize,
77 );
78
80 const G4Track& ,
81 const G4Step&
82 );
83
85 const G4Track&,
86 G4double previousStepSize,
87 G4double currentMinimumStep,
88 G4double& currentSafety,
89 G4GPILSelection* selection
90 ) ;
91
93 const G4Track& ,
94 const G4Step&
95 );
96
97 // no operation in AtRestDoIt
99 const G4Track& ,
101 ) { return -1.0; };
102
103 // no operation in AtRestDoIt
105 const G4Track& ,
106 const G4Step&
107 ) {return 0;};
108
109 protected: // with description
111 G4double previousStepSize,
112 G4double currentMinimumStep,
113 G4double& currentSafety
114 )=0;
115 private:
116 // this is the returnd value of G4GPILSelection in
117 // the arguments of AlongStepGPIL()
118 G4GPILSelection valueGPILSelection;
119
120 protected: // with description
121 //------------------------------------------------------
123 G4double previousStepSize) ;
124
125 // these two methods are set/get methods for valueGPILSelection
127 { valueGPILSelection = selection;};
128
129 G4GPILSelection GetGPILSelection() const{return valueGPILSelection;};
130
131 private:
132 // hide default constructor and assignment operator as private
135
136 protected:
140};
141
142// -----------------------------------------
143// inlined function members implementation
144// -----------------------------------------
145
148 G4double
149 )
150{
151 // dummy routine
152 ;
153}
154
155
157 const G4Track& ,
158 const G4Step&
159 )
160{
161// clear NumberOfInteractionLengthLeft
163 return pParticleChange;
164}
165
167 const G4Track& ,
168 const G4Step&
169 )
170{
171// clear NumberOfInteractionLengthLeft
173 return pParticleChange;
174}
175
177 const G4Track& track,
178 G4double previousStepSize,
179 G4double currentMinimumStep,
180 G4double& currentSafety,
181 G4GPILSelection* selection
182 )
183{
184 // GPILSelection is set to defaule value of CandidateForSelection
185 valueGPILSelection = CandidateForSelection;
186
187 // get Step limit proposed by the process
188 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
189
190 // set return value for G4GPILSelection
191 *selection = valueGPILSelection;
192
193 if (verboseLevel>1){
194 G4cout << "G4IVContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength ";
195 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
196 track.GetDynamicParticle()->DumpInfo();
197 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
198 G4cout << "IntractionLength= " << steplength/CLHEP::cm <<"[cm] " <<G4endl;
199 }
200 return steplength ;
201}
202
203#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4GPILSelection
@ CandidateForSelection
G4ProcessType
@ fNotDefined
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void DumpInfo(G4int mode=0) const
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
virtual G4double GetContinuousStepLimit(const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
void SetGPILSelection(G4GPILSelection selection)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
const G4String & GetName() const
Definition: G4Material.hh:177
Definition: G4Step.hh:78
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:418
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379