Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4WrapperProcess.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// G4WrapperProcess class implementation
27//
28// Author: H.Kurahige, 18 December 1996
29// --------------------------------------------------------------------
30
31#include "G4WrapperProcess.hh"
32
33// --------------------------------------------------------------------
35 G4ProcessType aType)
36 : G4VProcess(aName, aType)
37{
38}
39
40// --------------------------------------------------------------------
42 : G4VProcess(*((G4VProcess*)(&right))), pRegProcess(right.pRegProcess)
43{
44}
45
46// --------------------------------------------------------------------
48{
49}
50
51// --------------------------------------------------------------------
53{
55}
56
57// --------------------------------------------------------------------
60 G4double previousStepSize,
61 G4double currentMinimumStep,
62 G4double& proposedSafety,
63 G4GPILSelection* selection )
64{
65 return pRegProcess->
67 previousStepSize,
68 currentMinimumStep,
69 proposedSafety,
70 selection );
71}
72
73// --------------------------------------------------------------------
77{
79}
80
81// --------------------------------------------------------------------
84 G4double previousStepSize,
86{
88 previousStepSize,
89 condition );
90}
91
92// --------------------------------------------------------------------
94{
96}
97
99{
101}
102
103// --------------------------------------------------------------------
105 const G4Step& stepData )
106{
107 return pRegProcess->PostStepDoIt( track, stepData );
108}
109
110// --------------------------------------------------------------------
112 const G4Step& stepData )
113{
114 return pRegProcess->AlongStepDoIt( track, stepData );
115}
116
117// --------------------------------------------------------------------
119 const G4Step& stepData )
120{
121 return pRegProcess->AtRestDoIt( track, stepData );
122}
123
124// --------------------------------------------------------------------
126{
127 return pRegProcess->IsApplicable(particle);
128}
129
130// --------------------------------------------------------------------
132{
133 return pRegProcess->BuildPhysicsTable(particle);
134}
135
136// --------------------------------------------------------------------
138{
139 return pRegProcess->PreparePhysicsTable(particle);
140}
141
142// --------------------------------------------------------------------
145 const G4String& directory,
146 G4bool ascii)
147{
148 return pRegProcess->StorePhysicsTable(particle, directory, ascii);
149}
150
151// --------------------------------------------------------------------
154 const G4String& directory,
155 G4bool ascii)
156{
157 return pRegProcess->RetrievePhysicsTable(particle, directory, ascii);
158}
159
160// --------------------------------------------------------------------
162{
164}
165
166// --------------------------------------------------------------------
168{
170}
171
172// --------------------------------------------------------------------
174{
175 pRegProcess = process;
176 theProcessName += process->GetProcessName();
177 theProcessType = process->GetProcessType();
178}
179
180// --------------------------------------------------------------------
182{
183 return pRegProcess;
184}
185
186// --------------------------------------------------------------------
188{
189 G4WrapperProcess* master = static_cast<G4WrapperProcess*>(masterP);
190 // Cannot use getter because it returns "const" and we do not want
191 // to cast away constness explicitly (even if this is the same)
193}
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4GPILSelection
G4ProcessType
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
Definition: G4Step.hh:62
G4ProcessType theProcessType
Definition: G4VProcess.hh:350
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:80
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:211
virtual const G4ProcessManager * GetProcessManager()
Definition: G4VProcess.hh:498
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:392
virtual void SetMasterProcess(G4VProcess *masterP)
Definition: G4VProcess.cc:218
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:182
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:87
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:187
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:206
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:194
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)=0
G4String theProcessName
Definition: G4VProcess.hh:345
virtual void SetProcessManager(const G4ProcessManager *)
Definition: G4VProcess.hh:492
virtual void EndTracking()
Definition: G4VProcess.cc:102
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
virtual void SetMasterProcess(G4VProcess *masterP)
virtual void SetProcessManager(const G4ProcessManager *)
virtual void EndTracking()
virtual G4bool IsApplicable(const G4ParticleDefinition &)
virtual const G4VProcess * GetRegisteredProcess() const
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual void RegisterProcess(G4VProcess *)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
virtual ~G4WrapperProcess()
virtual void StartTracking(G4Track *)
virtual const G4ProcessManager * GetProcessManager()
G4VProcess * pRegProcess
virtual void ResetNumberOfInteractionLengthLeft()
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
G4WrapperProcess(const G4String &aName="Wrapped", G4ProcessType aType=fNotDefined)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)