Geant4  10.01.p02
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 85244 2014-10-27 08:24:13Z 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, kara@cenbg.in2p3.fr
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 
53 #ifndef G4VITRESTDISCRETEPROCESS_H
54 #define G4VITRESTDISCRETEPROCESS_H
55 
56 #include <CLHEP/Units/SystemOfUnits.h>
57 
58 #include "G4VITProcess.hh"
59 
65 {
66  // Abstract class which defines the public behavior of
67  // rest + discrete physics interactions.
68  public:
69 
71  G4ProcessType aType = fNotDefined );
73 
74  virtual ~G4VITRestDiscreteProcess();
75 
76  public :// with description
78  const G4Track& track,
79  G4double previousStepSize,
81  );
82 
84  const G4Track& ,
85  const G4Step&
86  );
87 
89  const G4Track& ,
91  );
92 
94  const G4Track& ,
95  const G4Step&
96  );
97 
98  // no operation in AlongStepDoIt
100  const G4Track&,
101  G4double ,
102  G4double ,
103  G4double& ,
105  ){ return -1.0; }
106 
107  // no operation in AlongStepDoIt
109  const G4Track& ,
110  const G4Step&
111  ) {return 0;}
112 
113  protected:// with description
114  virtual G4double GetMeanFreePath(const G4Track& aTrack,
115  G4double previousStepSize,
116  G4ForceCondition* condition
117  )=0;
118  // Calculates from the macroscopic cross section a mean
119  // free path, the value is returned in units of distance.
120 
121  virtual G4double GetMeanLifeTime(const G4Track& aTrack,G4ForceCondition* condition)=0;
122  // Calculates the mean life-time (i.e. for decays) of the
123  // particle at rest due to the occurence of the given process,
124  // or converts the probability of interaction (i.e. for
125  // annihilation) into the life-time of the particle for the
126  // occurence of the given process.
127 
128  private:
129  // hide default constructor and assignment operator as private
132 
133 };
134 
135 // -----------------------------------------
136 // inlined function members implementation
137 // -----------------------------------------
138 inline
140  const G4Track& track,
141  G4double previousStepSize,
143  )
144 {
145  if ( (previousStepSize < 0.0) || (fpState->theNumberOfInteractionLengthLeft<=0.0)) {
146  // beggining of tracking (or just after DoIt of this process)
148  } else if ( previousStepSize > 0.0) {
149  // subtract NumberOfInteractionLengthLeft
150  SubtractNumberOfInteractionLengthLeft(previousStepSize);
151  } else {
152  // zero step
153  // DO NOTHING
154  }
155 
156  // condition is set to "Not Forced"
157  *condition = NotForced;
158 
159  // get mean free path
160  fpState->currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
161 
162  G4double value;
163  if (fpState->currentInteractionLength <DBL_MAX) {
164  value = fpState->theNumberOfInteractionLengthLeft * (fpState->currentInteractionLength);
165  } else {
166  value = DBL_MAX;
167  }
168 #ifdef G4VERBOSE
169  if (verboseLevel>1){
170  G4cout << "G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength ";
171  G4cout << "[ " << GetProcessName() << "]" <<G4endl;
172  track.GetDynamicParticle()->DumpInfo();
173  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
174  G4cout << "InteractionLength= " << value/CLHEP::cm <<"[cm] " <<G4endl;
175  }
176 #endif
177  return value;
178 }
179 
181  const G4Track& ,
182  const G4Step&
183  )
184 {
185 // reset NumberOfInteractionLengthLeft
187 
188  return pParticleChange;
189 }
190 
192  const G4Track& track,
194  )
195 {
196  // beggining of tracking
198 
199  // condition is set to "Not Forced"
200  *condition = NotForced;
201 
202  // get mean life time
203  fpState->currentInteractionLength = GetMeanLifeTime(track, condition);
204 
205 #ifdef G4VERBOSE
206  if ((fpState->currentInteractionLength <0.0) || (verboseLevel>2)){
207  G4cout << "G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength ";
208  G4cout << "[ " << GetProcessName() << "]" <<G4endl;
209  track.GetDynamicParticle()->DumpInfo();
210  G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
211  G4cout << "MeanLifeTime = " << fpState->currentInteractionLength/CLHEP::ns << "[ns]" <<G4endl;
212  }
213 #endif
214 
215  return fpState->theNumberOfInteractionLengthLeft * (fpState->currentInteractionLength);
216 }
217 
218 
220  const G4Track&,
221  const G4Step&
222  )
223 {
224 // clear NumberOfInteractionLengthLeft
226 
227  return pParticleChange;
228 }
229 
230 
231 
232 #endif
233 
234 
static const double cm
Definition: G4SIunits.hh:106
G4double condition(const G4ErrorSymMatrix &m)
virtual void ClearNumberOfInteractionLengthLeft()
Identical to G4VRestDiscreteProcess with dependency from G4VITProcess.
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4DynamicParticle * GetDynamicParticle() const
virtual void ResetNumberOfInteractionLengthLeft()
WARNING : Redefine the method of G4VProcess reset (determine the value of)NumberOfInteractionLengthLe...
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
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
G4VITProcess inherits from G4VProcess.
Definition: G4VITProcess.hh:99
#define DBL_MAX
Definition: templates.hh:83
#define ns
Definition: xmlparse.cc:597
G4GPILSelection
G4ProcessType
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)