Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VRestContinuousDiscreteProcess.cc
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: G4VRestContinuousDiscreteProcess.cc 71231 2013-06-12 13:06:28Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------
31 // GEANT 4 class implementation file
32 //
33 // History: first implementation, based on object model of
34 // 2nd December 1995, G.Cosmo
35 // --------------------------------------------------------------
36 // New Physics scheme 8 Jan. 1997 H.Kurahige
37 // ------------------------------------------------------------
38 
40 #include "G4SystemOfUnits.hh"
41 
42 #include "G4Step.hh"
43 #include "G4Track.hh"
44 #include "G4MaterialTable.hh"
45 #include "G4VParticleChange.hh"
46 
48  :G4VProcess("No Name Discrete Process"),
49  valueGPILSelection(CandidateForSelection)
50 {
51  G4Exception("G4VRestContinuousDiscreteProcess::G4VRestContinuousDiscreteProcess()",
52  "ProcMan102",JustWarning,"Default constructor is called");
53 }
54 
56  : G4VProcess(aName, aType),
57  valueGPILSelection(CandidateForSelection)
58 {
59 }
60 
62 {
63 }
64 
66  : G4VProcess(right),
67  valueGPILSelection(right.valueGPILSelection)
68 {
69 }
70 
72  const G4Track& track,
74  )
75 {
76  // beggining of tracking
78 
79  // condition is set to "Not Forced"
80  *condition = NotForced;
81 
82  // get mean life time
83  currentInteractionLength = GetMeanLifeTime(track, condition);
84 
85 #ifdef G4VERBOSE
86  if ((currentInteractionLength <0.0) || (verboseLevel>2)){
87  G4cout << "G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength ";
88  G4cout << "[ " << GetProcessName() << "]" <<G4endl;
89  track.GetDynamicParticle()->DumpInfo();
90  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
91  G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl;
92  }
93 #endif
94 
96 }
97 
98 
100  const G4Track&,
101  const G4Step&
102  )
103 {
104 // clear NumberOfInteractionLengthLeft
106 
107  return pParticleChange;
108 }
109 
111  const G4Track& track,
112  G4double previousStepSize,
113  G4double currentMinimumStep,
114  G4double& currentSafety,
115  G4GPILSelection* selection
116  )
117 {
118  // GPILSelection is set to defaule value of CandidateForSelection
119  valueGPILSelection = CandidateForSelection;
120 
121  // get Step limit proposed by the process
122  G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
123 
124  // set return value for G4GPILSelection
125  *selection = valueGPILSelection;
126 
127 #ifdef G4VERBOSE
128  if (verboseLevel>1){
129  G4cout << "G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength ";
130  G4cout << "[ " << GetProcessName() << "]" <<G4endl;
131  track.GetDynamicParticle()->DumpInfo();
132  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
133  G4cout << "IntractionLength= " << steplength/cm <<"[cm] " <<G4endl;
134  }
135 #endif
136  return steplength ;
137 }
138 
140  const G4Track& ,
141  const G4Step&
142  )
143 {
144  return pParticleChange;
145 }
146 
148  const G4Track& track,
149  G4double previousStepSize,
151  )
152 {
153  if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0)) {
154  // beggining of tracking (or just after DoIt of this process)
156  } else if ( previousStepSize > 0.0) {
157  // subtract NumberOfInteractionLengthLeft
158  SubtractNumberOfInteractionLengthLeft(previousStepSize);
159  } else {
160  // zero step
161  // DO NOTHING
162  }
163 
164  // condition is set to "Not Forced"
165  *condition = NotForced;
166 
167  // get mean free path
168  currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
169 
170 
171  G4double value;
174  } else {
175  value = DBL_MAX;
176  }
177 #ifdef G4VERBOSE
178  if (verboseLevel>1){
179  G4cout << "G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength ";
180  G4cout << "[ " << GetProcessName() << "]" <<G4endl;
181  track.GetDynamicParticle()->DumpInfo();
182  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
183  G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
184  }
185 #endif
186  return value;
187 }
188 
190  const G4Track& ,
191  const G4Step&
192  )
193 {
194 // clear NumberOfInteractionLengthLeft
196 
197  return pParticleChange;
198 }
G4double condition(const G4ErrorSymMatrix &m)
G4VRestContinuousDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4DynamicParticle * GetDynamicParticle() const
const G4String & GetName() const
Definition: G4Material.hh:178
void DumpInfo(G4int mode=0) const
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual G4double GetContinuousStepLimit(const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
static constexpr double cm
Definition: G4SIunits.hh:119
G4double currentInteractionLength
Definition: G4VProcess.hh:297
Definition: G4Step.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4Material * GetMaterial() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
#define G4endl
Definition: G4ios.hh:61
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
Definition: G4VProcess.hh:544
double G4double
Definition: G4Types.hh:76
G4ForceCondition
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
#define DBL_MAX
Definition: templates.hh:83
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
#define ns
Definition: xmlparse.cc:614
G4GPILSelection
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)=0
G4ProcessType