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
G4AdjointProcessEquivalentToDirectProcess.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// ------------------------------------------------------------
27//
28// This class is for adjoint process equivalent to direct process
29//
30// ------------------------------------------------------------
31// Created by L.Desorgher 25 Sept. 2009
32// ------------------------------------------------------------
33
35
36#include "G4DynamicParticle.hh"
38#include "G4VProcess.hh"
39
42 const G4String& aName, G4VProcess* aProcess,
43 G4ParticleDefinition* fwd_particle_def)
44 : G4VProcess(aName)
45{
46 fDirectProcess = aProcess;
47 theProcessType = fDirectProcess->GetProcessType();
48 fFwdParticleDef = fwd_particle_def;
49}
50
53{
54 if(fDirectProcess != nullptr)
55 delete fDirectProcess;
56}
57
60{
61 fDirectProcess->ResetNumberOfInteractionLengthLeft();
62}
63
66 G4double previousStepSize,
67 G4double currentMinimumStep,
68 G4double& proposedSafety,
69 G4GPILSelection* selection)
70{
71 // Change the particle definition to the direct one
72 G4DynamicParticle* theDynPart =
73 const_cast<G4DynamicParticle*>(track.GetDynamicParticle());
74 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
75
76 G4DecayProducts* decayProducts =
77 const_cast<G4DecayProducts*>(theDynPart->GetPreAssignedDecayProducts());
79 theDynPart->SetDefinition(fFwdParticleDef);
80
81 // Call the direct process
83 track, previousStepSize, currentMinimumStep, proposedSafety, selection);
84
85 // Restore the adjoint particle definition to the direct one
86 theDynPart->SetDefinition(adjPartDef);
87 theDynPart->SetPreAssignedDecayProducts(decayProducts);
88
89 return GPIL;
90}
91
94 const G4Track& track, G4ForceCondition* condition)
95{
96 // Change the particle definition to the direct one
97 G4DynamicParticle* theDynPart =
98 const_cast<G4DynamicParticle*>(track.GetDynamicParticle());
99 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
100
101 G4DecayProducts* decayProducts =
102 const_cast<G4DecayProducts*>(theDynPart->GetPreAssignedDecayProducts());
104 theDynPart->SetDefinition(fFwdParticleDef);
105
106 // Call the direct process
107 G4double GPIL =
108 fDirectProcess->AtRestGetPhysicalInteractionLength(track, condition);
109
110 // Restore the adjoint particle definition to the direct one
111 theDynPart->SetDefinition(adjPartDef);
112 theDynPart->SetPreAssignedDecayProducts(decayProducts);
113
114 return GPIL;
115}
116
119 const G4Track& track, G4double previousStepSize, G4ForceCondition* condition)
120{
121 // Change the particle definition to the direct one
122 G4DynamicParticle* theDynPart =
123 const_cast<G4DynamicParticle*>(track.GetDynamicParticle());
124 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
125
126 G4DecayProducts* decayProducts =
127 const_cast<G4DecayProducts*>(theDynPart->GetPreAssignedDecayProducts());
128
130 theDynPart->SetDefinition(fFwdParticleDef);
131
132 // Call the direct process
133 G4double GPIL = fDirectProcess->PostStepGetPhysicalInteractionLength(
134 track, previousStepSize, condition);
135
136 // Restore the adjoint particle definition to the direct one
137 theDynPart->SetDefinition(adjPartDef);
138 theDynPart->SetPreAssignedDecayProducts(decayProducts);
139
140 return GPIL;
141}
142
144 const G4Track& track, const G4Step& stepData)
145{
146 // Change the particle definition to the direct one
147 G4DynamicParticle* theDynPart =
148 const_cast<G4DynamicParticle*>(track.GetDynamicParticle());
149 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
150
151 G4DecayProducts* decayProducts =
152 const_cast<G4DecayProducts*>(theDynPart->GetPreAssignedDecayProducts());
153
155 theDynPart->SetDefinition(fFwdParticleDef);
156
157 // Call the direct process
158 G4VParticleChange* partChange = fDirectProcess->PostStepDoIt(track, stepData);
159
160 // Restore the adjoint particle definition to the direct one
161 theDynPart->SetDefinition(adjPartDef);
162 theDynPart->SetPreAssignedDecayProducts(decayProducts);
163
164 return partChange;
165}
166
168 const G4Track& track, const G4Step& stepData)
169{
170 // Change the particle definition to the direct one
171 G4DynamicParticle* theDynPart =
172 const_cast<G4DynamicParticle*>(track.GetDynamicParticle());
173 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
174
175 G4DecayProducts* decayProducts =
176 const_cast<G4DecayProducts*>(theDynPart->GetPreAssignedDecayProducts());
177
179 theDynPart->SetDefinition(fFwdParticleDef);
180
181 // Call the direct process
182 G4VParticleChange* partChange =
183 fDirectProcess->AlongStepDoIt(track, stepData);
184
185 // Restore the adjoint particle definition to the direct one
186 theDynPart->SetDefinition(adjPartDef);
187 theDynPart->SetPreAssignedDecayProducts(decayProducts);
188
189 return partChange;
190}
191
193 const G4Track& track, const G4Step& stepData)
194{
195 // Change the particle definition to the direct one
196 G4DynamicParticle* theDynPart =
197 const_cast<G4DynamicParticle*>(track.GetDynamicParticle());
198 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
199
200 G4DecayProducts* decayProducts =
201 const_cast<G4DecayProducts*>(theDynPart->GetPreAssignedDecayProducts());
202
204 theDynPart->SetDefinition(fFwdParticleDef);
205
206 // Call the direct process
207 G4VParticleChange* partChange = fDirectProcess->AtRestDoIt(track, stepData);
208
209 // Restore the adjoint particle definition to the direct one
210 theDynPart->SetDefinition(adjPartDef);
211 theDynPart->SetPreAssignedDecayProducts(decayProducts);
212
213 return partChange;
214}
215
218{
219 return fDirectProcess->IsApplicable(*fFwdParticleDef);
220}
221
224{
225 return fDirectProcess->BuildPhysicsTable(*fFwdParticleDef);
226}
227
230{
231 return fDirectProcess->PreparePhysicsTable(*fFwdParticleDef);
232}
233
235 const G4ParticleDefinition*, const G4String& directory, G4bool ascii)
236{
237 return fDirectProcess->StorePhysicsTable(fFwdParticleDef, directory, ascii);
238}
239
241 const G4ParticleDefinition*, const G4String& directory, G4bool ascii)
242{
243 return fDirectProcess->RetrievePhysicsTable(fFwdParticleDef, directory,
244 ascii);
245}
246
248{
249 // Change the particle definition to the direct one
250 G4DynamicParticle* theDynPart =
251 const_cast<G4DynamicParticle*>(track->GetDynamicParticle());
252 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
253
254 G4DecayProducts* decayProducts =
255 const_cast<G4DecayProducts*>(theDynPart->GetPreAssignedDecayProducts());
257 theDynPart->SetDefinition(fFwdParticleDef);
258
259 fDirectProcess->StartTracking(track);
260
261 // Restore the adjoint particle definition to the direct one
262 theDynPart->SetDefinition(adjPartDef);
263 theDynPart->SetPreAssignedDecayProducts(decayProducts);
264
265 return;
266}
267
269{
270 fDirectProcess->EndTracking();
271}
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4GPILSelection
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData) override
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData) override
void PreparePhysicsTable(const G4ParticleDefinition &) override
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition) override
G4AdjointProcessEquivalentToDirectProcess(const G4String &aName, G4VProcess *aProcess, G4ParticleDefinition *fwd_particle_def)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData) override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4bool IsApplicable(const G4ParticleDefinition &) override
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection) override
void SetPreAssignedDecayProducts(G4DecayProducts *aDecayProducts)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
const G4DecayProducts * GetPreAssignedDecayProducts() const
G4ParticleDefinition * GetDefinition() const
Definition: G4Step.hh:62
const G4DynamicParticle * GetDynamicParticle() const
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
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:392
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 G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=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
virtual void EndTracking()
Definition: G4VProcess.cc:102