Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 100802 2016-11-02 14:55:27Z gcosmo $
27 //
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
34 
35 // The code is developed in the framework of the ESA AO7146
36 //
37 // We would be very happy hearing from you, send us your feedback! :)
38 //
39 // In order for Geant4-DNA to be maintained and still open-source,
40 // article citations are crucial.
41 // If you use Geant4-DNA chemistry and you publish papers about your software,
42 // in addition to the general paper on Geant4-DNA:
43 //
44 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
45 //
46 // we would be very happy if you could please also cite the following
47 // reference papers on chemistry:
48 //
49 // J. Comput. Phys. 274 (2014) 841-882
50 // Prog. Nucl. Sci. Tec. 2 (2011) 503-508
51 
52 #ifndef G4VITRESTDISCRETEPROCESS_H
53 #define G4VITRESTDISCRETEPROCESS_H
54 
56 
57 #include "G4VITProcess.hh"
58 
64 {
65  // Abstract class which defines the public behavior of
66  // rest + discrete physics interactions.
67 public:
68 
71 
72  virtual ~G4VITRestDiscreteProcess();
73 
74 public:
75  // with description
76  virtual G4double
78  G4double previousStepSize,
80 
81  virtual G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
82 
85 
86  virtual G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&);
87 
88  // no operation in AlongStepDoIt
90  G4double,
91  G4double,
92  G4double&,
94  {
95  return -1.0;
96  }
97 
98  // no operation in AlongStepDoIt
99  virtual G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&)
100  {
101  return 0;
102  }
103 
104 protected:
105  // with description
106  virtual G4double GetMeanFreePath(const G4Track& aTrack,
107  G4double previousStepSize,
108  G4ForceCondition* condition)=0;
109  // Calculates from the macroscopic cross section a mean
110  // free path, the value is returned in units of distance.
111 
112  virtual G4double GetMeanLifeTime(const G4Track& aTrack,
113  G4ForceCondition* condition)=0;
114  // Calculates the mean life-time (i.e. for decays) of the
115  // particle at rest due to the occurence of the given process,
116  // or converts the probability of interaction (i.e. for
117  // annihilation) into the life-time of the particle for the
118  // occurence of the given process.
119 
120 protected:
121  // hide default constructor and assignment operator as private
124 
125 };
126 
127 // -----------------------------------------
128 // inlined function members implementation
129 // -----------------------------------------
131  G4double previousStepSize,
133 {
134  if((previousStepSize < 0.0) || (fpState->theNumberOfInteractionLengthLeft
135  <= 0.0))
136  {
137  // beggining of tracking (or just after DoIt of this process)
139  }
140  else if(previousStepSize > 0.0)
141  {
142  // subtract NumberOfInteractionLengthLeft
143  SubtractNumberOfInteractionLengthLeft(previousStepSize);
144  }
145  else
146  {
147  // zero step
148  // DO NOTHING
149  }
150 
151  // condition is set to "Not Forced"
152  *condition = NotForced;
153 
154  // get mean free path
155  fpState->currentInteractionLength = GetMeanFreePath(track,
156  previousStepSize,
157  condition);
158 
159  G4double value;
160  if(fpState->currentInteractionLength < DBL_MAX)
161  {
162  value = fpState->theNumberOfInteractionLengthLeft * (fpState->currentInteractionLength);
163  }
164  else
165  {
166  value = DBL_MAX;
167  }
168 #ifdef G4VERBOSE
169  if(verboseLevel > 1)
170  {
171  G4cout << "G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength ";
172  G4cout << "[ " << GetProcessName() << "]" << G4endl;
173  track.GetDynamicParticle()->DumpInfo();
174  G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
175  G4cout << "InteractionLength= " << value / CLHEP::cm << "[cm] " << G4endl;
176  }
177 #endif
178  return value;
179 }
180 
182  const G4Step&)
183 {
184 // reset NumberOfInteractionLengthLeft
186 
187  return pParticleChange;
188 }
189 
192 {
193  // beggining of tracking
195 
196  // condition is set to "Not Forced"
197  *condition = NotForced;
198 
199  // get mean life time
200  fpState->currentInteractionLength = GetMeanLifeTime(track, condition);
201 
202 #ifdef G4VERBOSE
203  if((fpState->currentInteractionLength < 0.0) || (verboseLevel > 2))
204  {
205  G4cout << "G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength ";
206  G4cout << "[ " << GetProcessName() << "]" << G4endl;
207  track.GetDynamicParticle()->DumpInfo();
208  G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
209  G4cout << "MeanLifeTime = " << fpState->currentInteractionLength / CLHEP::ns
210  << "[ns]" << G4endl;
211  }
212 #endif
213 
214  return fpState->theNumberOfInteractionLengthLeft
215  * (fpState->currentInteractionLength);
216 }
217 
219  const G4Step&)
220 {
221 // clear NumberOfInteractionLengthLeft
223 
224  return pParticleChange;
225 }
226 
227 #endif
228 
G4double condition(const G4ErrorSymMatrix &m)
virtual void ClearNumberOfInteractionLengthLeft()
Identical to G4VRestDiscreteProcess with dependency from G4VITProcess.
static constexpr double cm
Definition: SystemOfUnits.h:99
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4DynamicParticle * GetDynamicParticle() const
static constexpr double ns
virtual void ResetNumberOfInteractionLengthLeft()
const G4String & GetName() const
Definition: G4Material.hh:178
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void DumpInfo(G4int mode=0) const
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4VITRestDiscreteProcess & operator=(const G4VITRestDiscreteProcess &right)
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4shared_ptr< G4ProcessState > fpState
Definition: G4Step.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4Material * GetMaterial() const
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)=0
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
#define G4endl
Definition: G4ios.hh:61
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
double G4double
Definition: G4Types.hh:76
G4ForceCondition
#define DBL_MAX
Definition: templates.hh:83
G4GPILSelection
G4ProcessType
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)