Geant4  10.00.p03
G4VContinuousDiscreteProcess.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: G4VContinuousDiscreteProcess.hh 71231 2013-06-12 13:06:28Z gcosmo $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 // Abstract class which defines the public behavior of
34 //
35 // Class Description
36 // discrete physics interactions.
37 //
38 // ------------------------------------------------------------
39 // New Physics scheme 8 Mar. 1997 H.Kurahige
40 // ------------------------------------------------------------
41 // fix bugs in GetGPILSelection() 24 Jan. 1998 H.Kurashige
42 // modified for new ParticleChange 12 Mar. 1998 H.Kurashige
43 // Fixed a bug in PostStepGetPhysicalInteractionLength
44 // 15 Apr. 2002 H.Kurashige
45 //
46 
47 #ifndef G4VContinuousDiscreteProcess_h
48 #define G4VContinuousDiscreteProcess_h 1
49 
50 #include "globals.hh"
51 #include "G4ios.hh"
52 
53 #include "G4VProcess.hh"
54 
56 {
57  // Abstract class which defines the public behavior of
58  // discrete physics interactions.
59  public:
60 
62  G4ProcessType aType = fNotDefined );
64 
66 
67  public :// with description
69  const G4Track& track,
70  G4double previousStepSize,
72  );
73 
75  const G4Track& ,
76  const G4Step&
77  );
78 
80  const G4Track&,
81  G4double previousStepSize,
82  G4double currentMinimumStep,
83  G4double& currentSafety,
84  G4GPILSelection* selection
85  );
86 
88  const G4Track& ,
89  const G4Step&
90  );
91 
92  // no operation in AtRestDoIt
94  const G4Track& ,
96  ) { return -1.0; };
97 
98  // no operation in AtRestDoIt
100  const G4Track& ,
101  const G4Step&
102  ) {return 0;};
103 
104  protected:// with description
105  virtual G4double GetMeanFreePath(const G4Track& aTrack,
106  G4double previousStepSize,
107  G4ForceCondition* condition
108  )=0;
109  // Calculates from the macroscopic cross section a mean
110  // free path, the value is returned in units of distance.
111 
112  protected:// with description
113  virtual G4double GetContinuousStepLimit(const G4Track& aTrack,
114  G4double previousStepSize,
115  G4double currentMinimumStep,
116  G4double& currentSafety
117  )=0;
118  private:
119  // this is the returnd value of G4GPILSelection in
120  // the arguments of AlongStepGPIL()
122 
123  protected:// with description
124  // these two methods are set/get methods for valueGPILSelection
126  { valueGPILSelection = selection;};
127 
129 
130  private:
131  // hide default constructor and assignment operator as private
134 
135 };
136 
137 #endif
138 
139 
140 
141 
142 
143 
G4double condition(const G4ErrorSymMatrix &m)
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
void SetGPILSelection(G4GPILSelection selection)
G4VContinuousDiscreteProcess & operator=(const G4VContinuousDiscreteProcess &right)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
virtual G4double GetContinuousStepLimit(const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
Definition: G4Step.hh:76
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
double G4double
Definition: G4Types.hh:76
G4ForceCondition
G4GPILSelection
G4ProcessType