source: trunk/source/processes/management/include/G4VContinuousProcess.hh@ 831

Last change on this file since 831 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 6.7 KB
Line 
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: G4VContinuousProcess.hh,v 1.6 2006/06/29 21:07:46 gunter Exp $
28// GEANT4 tag $Name: $
29//
30//
31// ------------------------------------------------------------
32// GEANT 4 class header file
33//
34// History: first implementation, based on object model of
35// 2nd December 1995, G.Cosmo
36// add G4VContinuousProcess(const G4String&) 24 Jul 1996, Hisaya kurashige
37//
38// Class Description
39// Abstract class which defines the public behavior of
40// Continuous physics interactions.
41//
42// ------------------------------------------------------------
43// New Physics scheme 18 Dec. 1996 H.Kurahige
44// ------------------------------------------------------------
45// modified 25 Feb. 1997 H.Kurahige
46// modified 8 Mar. 1997 H.Kurashige
47// modified 22 Mar. 1997 H.Kurashige
48// modified 26 Mar. 1997 H.Kurashige
49// modified AlongStepGPIL etc. 17 Dec. 1997 H.Kurashige
50// fix bugs in GetGPILSelection() 24 Jan. 1998 H.Kurashige
51// modified for new ParticleChange 12 Mar. 1998 H.Kurashige
52
53#ifndef G4VContinuousProcess_h
54#define G4VContinuousProcess_h 1
55
56#include "globals.hh"
57#include "G4ios.hh"
58
59#include "G4VProcess.hh"
60
61class G4VContinuousProcess : public G4VProcess
62{
63 // Abstract class which defines the public behavior of
64 // Continuous physics interactions.
65 public:
66
67 G4VContinuousProcess(const G4String& ,
68 G4ProcessType aType = fNotDefined );
69 G4VContinuousProcess(G4VContinuousProcess &);
70
71 virtual ~G4VContinuousProcess();
72
73 public: // with description
74 virtual G4double AlongStepGetPhysicalInteractionLength(
75 const G4Track& track,
76 G4double previousStepSize,
77 G4double currentMinimumStep,
78 G4double& proposedSafety,
79 G4GPILSelection* selection
80 );
81
82 virtual G4VParticleChange* AlongStepDoIt(
83 const G4Track& ,
84 const G4Step&
85 );
86
87 // no operation in AtRestDoIt and PostStepDoIt
88 virtual G4double PostStepGetPhysicalInteractionLength(
89 const G4Track&,
90 G4double,
91 G4ForceCondition*
92 ){ return -1.0; };
93
94 virtual G4double AtRestGetPhysicalInteractionLength(
95 const G4Track& ,
96 G4ForceCondition*
97 ) { return -1.0; };
98
99 // no operation in AtRestDoIt and PostStepDoIt
100 virtual G4VParticleChange* AtRestDoIt(
101 const G4Track& ,
102 const G4Step&
103 ) {return 0;};
104
105 virtual G4VParticleChange* PostStepDoIt(
106 const G4Track& ,
107 const G4Step&
108 ) {return 0;};
109
110 protected: // with description
111 virtual G4double GetContinuousStepLimit(const G4Track& aTrack,
112 G4double previousStepSize,
113 G4double currentMinimumStep,
114 G4double& currentSafety
115 )=0;
116 // This pure virtual function is used to calculate step limit
117 // for AlongStep in the derived processes
118
119 private:
120 // this is the returnd value of G4GPILSelection in
121 // the arguments of AlongStepGPIL()
122 G4GPILSelection valueGPILSelection;
123
124 protected: // with description
125 // these two methods are set/get methods for valueGPILSelection
126 void SetGPILSelection(G4GPILSelection selection)
127 { valueGPILSelection = selection;};
128
129 G4GPILSelection GetGPILSelection() const{return valueGPILSelection;};
130
131
132 private:
133 // hide default constructor and assignment operator as private
134 G4VContinuousProcess();
135 G4VContinuousProcess & operator=(const G4VContinuousProcess &right);
136
137};
138// -----------------------------------------
139// inlined function members implementation
140// -----------------------------------------
141#include "G4Step.hh"
142#include "G4Track.hh"
143#include "G4MaterialTable.hh"
144#include "G4VParticleChange.hh"
145
146inline G4double G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(
147 const G4Track& track,
148 G4double previousStepSize,
149 G4double currentMinimumStep,
150 G4double& currentSafety,
151 G4GPILSelection* selection
152 )
153{
154 // GPILSelection is set to defaule value of CandidateForSelection
155 valueGPILSelection = CandidateForSelection;
156
157 // get Step limit proposed by the process
158 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
159
160 // set return value for G4GPILSelection
161 *selection = valueGPILSelection;
162
163#ifdef G4VERBOSE
164 if (verboseLevel>1){
165 G4cout << "G4VContinuousProcess::AlongStepGetPhysicalInteractionLength ";
166 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
167 track.GetDynamicParticle()->DumpInfo();
168 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
169 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " <<G4endl;
170 }
171#endif
172
173 return steplength ;
174}
175
176inline G4VParticleChange* G4VContinuousProcess::AlongStepDoIt(
177 const G4Track& ,
178 const G4Step&
179 )
180{
181 return pParticleChange;
182}
183
184#endif
185
186
187
188
Note: See TracBrowser for help on using the repository browser.