Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VProcess.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// G4VProcess class implementation
27//
28// Authors:
29// - 2 December 1995, G.Cosmo - First implementation, based on object model
30// - 18 December 1996, H.Kurashige - New Physics scheme
31// --------------------------------------------------------------------
32
33#include "G4VProcess.hh"
34
36#include "G4SystemOfUnits.hh"
37
38#include "G4ProcessTable.hh"
39#include "G4PhysicsTable.hh"
40#include "G4MaterialTable.hh"
41#include "G4ElementTable.hh"
42#include "G4ElementVector.hh"
43#include "G4Log.hh"
44
45// --------------------------------------------------------------------
47 : theProcessName(aName), theProcessType(aType)
48{
50 fProcessTable = G4ProcessTable::GetProcessTable();
51 fProcessTable->RegisterProcess(this);
52}
53
54// --------------------------------------------------------------------
55G4VProcess::G4VProcess()
56{
57}
58
59// --------------------------------------------------------------------
61{
62 fProcessTable->DeRegisterProcess(this);
63}
64
65// --------------------------------------------------------------------
67 : theProcessName(right.theProcessName),
68 theProcessType(right.theProcessType),
69 theProcessSubType(right.theProcessSubType),
70 verboseLevel(right.verboseLevel),
71 enableAtRestDoIt(right.enableAtRestDoIt),
72 enableAlongStepDoIt(right.enableAlongStepDoIt),
73 enablePostStepDoIt(right.enablePostStepDoIt),
74 masterProcessShadow(right.masterProcessShadow),
75 fProcessTable(right.fProcessTable)
76{
77}
78
79// --------------------------------------------------------------------
81{
84}
85
86// --------------------------------------------------------------------
88{
92#ifdef G4VERBOSE
93 if (verboseLevel>2)
94 {
95 G4cout << "G4VProcess::StartTracking() - [" << theProcessName << "]"
96 << G4endl;
97 }
98#endif
99}
100
101// --------------------------------------------------------------------
103{
104#ifdef G4VERBOSE
105 if (verboseLevel>2)
106 {
107 G4cout << "G4VProcess::EndTracking() - [" << theProcessName << "]"
108 << G4endl;
109 }
110#endif
114}
115
116// --------------------------------------------------------------------
117namespace
118{
119 static const G4String typeNotDefined = "NotDefined";
120 static const G4String typeTransportation = "Transportation";
121 static const G4String typeElectromagnetic = "Electromagnetic";
122 static const G4String typeOptical = "Optical";
123 static const G4String typeHadronic = "Hadronic";
124 static const G4String typePhotolepton_hadron = "Photolepton_hadron";
125 static const G4String typeDecay = "Decay";
126 static const G4String typeGeneral = "General";
127 static const G4String typeParameterisation = "Parameterisation";
128 static const G4String typeUserDefined = "UserDefined";
129 static const G4String typePhonon = "Phonon";
130 static const G4String noType = "------";
131}
132
133// --------------------------------------------------------------------
135{
136 switch (aType)
137 {
138 case fNotDefined: return typeNotDefined; break;
139 case fTransportation: return typeTransportation; break;
140 case fElectromagnetic: return typeElectromagnetic; break;
141 case fOptical: return typeOptical; break;
142 case fHadronic: return typeHadronic; break;
143 case fPhotolepton_hadron: return typePhotolepton_hadron; break;
144 case fDecay: return typeDecay; break;
145 case fGeneral: return typeGeneral; break;
146 case fParameterisation: return typeParameterisation; break;
147 case fUserDefined: return typeUserDefined; break;
148 case fPhonon: return typePhonon; break;
149 default: ;
150 }
151 return noType;
152}
153
154// --------------------------------------------------------------------
156{
157 return this;
158}
159
160// --------------------------------------------------------------------
162{
163 return (this == &right);
164}
165
166// --------------------------------------------------------------------
168{
169 return (this != &right);
170}
171
172// --------------------------------------------------------------------
174{
175 G4cout << "Process Name " << theProcessName ;
176 G4cout << " : Type[" << GetProcessTypeName(theProcessType) << "]";
177 G4cout << " : SubType[" << theProcessSubType << "]"<< G4endl;
178}
179
180// --------------------------------------------------------------------
181void G4VProcess::ProcessDescription(std::ostream& outFile) const
182{
183 outFile << "This process has not yet been described\n";
184}
185
186// --------------------------------------------------------------------
189 const G4String& directory,
190 const G4String& tableName,
191 G4bool ascii )
192{
193 G4String thePhysicsTableFileExt;
194 if (ascii) thePhysicsTableFileExt = ".asc";
195 else thePhysicsTableFileExt = ".dat";
196
197 thePhysicsTableFileName = directory + "/";
198 thePhysicsTableFileName += tableName + "." + theProcessName + ".";
200 + thePhysicsTableFileExt;
201
203}
204
205// --------------------------------------------------------------------
207{
208 BuildPhysicsTable(part);
209}
210
211// --------------------------------------------------------------------
213{
215}
216
217// --------------------------------------------------------------------
219{
220 masterProcessShadow = masterP;
221}
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4ProcessType
@ fOptical
@ fPhonon
@ fParameterisation
@ fGeneral
@ fDecay
@ fElectromagnetic
@ fHadronic
@ fUserDefined
@ fTransportation
@ fPhotolepton_hadron
@ fNotDefined
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
const G4String & GetParticleName() const
static G4ProcessTable * GetProcessTable()
void RegisterProcess(G4VProcess *)
void DeRegisterProcess(G4VProcess *)
G4double currentInteractionLength
Definition: G4VProcess.hh:339
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:342
G4bool operator==(const G4VProcess &right) const
Definition: G4VProcess.cc:161
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:134
G4ProcessType theProcessType
Definition: G4VProcess.hh:350
virtual void ProcessDescription(std::ostream &outfile) const
Definition: G4VProcess.cc:181
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:80
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:46
virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.cc:212
G4String thePhysicsTableFileName
Definition: G4VProcess.hh:348
virtual void SetMasterProcess(G4VProcess *masterP)
Definition: G4VProcess.cc:218
virtual ~G4VProcess()
Definition: G4VProcess.cc:60
G4int verboseLevel
Definition: G4VProcess.hh:360
G4bool operator!=(const G4VProcess &right) const
Definition: G4VProcess.cc:167
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:335
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:87
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:187
G4int theProcessSubType
Definition: G4VProcess.hh:353
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:194
virtual const G4VProcess * GetCreatorProcess() const
Definition: G4VProcess.cc:155
virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition &part)
Definition: G4VProcess.cc:206
G4String theProcessName
Definition: G4VProcess.hh:345
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
virtual void DumpInfo() const
Definition: G4VProcess.cc:173
virtual void EndTracking()
Definition: G4VProcess.cc:102
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:188