Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointProcessEquivalentToDirectProcess Class Reference

#include <G4AdjointProcessEquivalentToDirectProcess.hh>

+ Inheritance diagram for G4AdjointProcessEquivalentToDirectProcess:

Public Member Functions

 G4AdjointProcessEquivalentToDirectProcess (const G4String &aName, G4VProcess *aProcess, G4ParticleDefinition *fwd_particle_def)
 
virtual ~G4AdjointProcessEquivalentToDirectProcess ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 55 of file G4AdjointProcessEquivalentToDirectProcess.hh.

Constructor & Destructor Documentation

◆ G4AdjointProcessEquivalentToDirectProcess()

G4AdjointProcessEquivalentToDirectProcess::G4AdjointProcessEquivalentToDirectProcess ( const G4String aName,
G4VProcess aProcess,
G4ParticleDefinition fwd_particle_def 
)

Definition at line 43 of file G4AdjointProcessEquivalentToDirectProcess.cc.

46:G4VProcess(aName)
47{
48 theDirectProcess =aProcess;
49 theProcessType = theDirectProcess->GetProcessType();
50 theFwdParticleDef = fwd_particle_def;
51}
G4ProcessType theProcessType
Definition: G4VProcess.hh:340
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:385

◆ ~G4AdjointProcessEquivalentToDirectProcess()

G4AdjointProcessEquivalentToDirectProcess::~G4AdjointProcessEquivalentToDirectProcess ( )
virtual

Definition at line 54 of file G4AdjointProcessEquivalentToDirectProcess.cc.

55{
56 if (theDirectProcess!=0) delete theDirectProcess;
57}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4AdjointProcessEquivalentToDirectProcess::AlongStepDoIt ( const G4Track track,
const G4Step stepData 
)
virtual

Implements G4VProcess.

Definition at line 208 of file G4AdjointProcessEquivalentToDirectProcess.cc.

210{
211 //Change the particle definition to the direct one
212 //------------------------------------------------
213 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
214 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
215
216 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
217
219 theDynPart->SetDefinition(theFwdParticleDef);
220
221
222 //Call the direct process
223 //----------------------
224 G4VParticleChange* partChange =theDirectProcess->AlongStepDoIt( track, stepData );
225
226 //Restore the adjoint particle definition to the direct one
227 //------------------------------------------------
228 theDynPart->SetDefinition(adjPartDef);
229 theDynPart->SetPreAssignedDecayProducts(decayProducts);
230
231 return partChange;
232}
void SetPreAssignedDecayProducts(G4DecayProducts *aDecayProducts)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
const G4DecayProducts * GetPreAssignedDecayProducts() const
G4ParticleDefinition * GetDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0

◆ AlongStepGetPhysicalInteractionLength()

G4double G4AdjointProcessEquivalentToDirectProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 64 of file G4AdjointProcessEquivalentToDirectProcess.cc.

70{
71
72
73 //Change the particle definition to the direct one
74 //------------------------------------------------
75 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
76 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
77
78 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
80 theDynPart->SetDefinition(theFwdParticleDef);
81
82
83 //Call the direct process
84 //----------------------
85 G4double GPIL = theDirectProcess->
87 previousStepSize,
88 currentMinimumStep,
89 proposedSafety,
90 selection );
91
92
93 //Restore the adjoint particle definition to the direct one
94 //------------------------------------------------
95 theDynPart->SetDefinition(adjPartDef);
96 theDynPart->SetPreAssignedDecayProducts(decayProducts);
97
98
99 return GPIL;
100
101}
double G4double
Definition: G4Types.hh:64
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)

Referenced by AlongStepGetPhysicalInteractionLength().

◆ AtRestDoIt()

G4VParticleChange * G4AdjointProcessEquivalentToDirectProcess::AtRestDoIt ( const G4Track track,
const G4Step stepData 
)
virtual

Implements G4VProcess.

Definition at line 234 of file G4AdjointProcessEquivalentToDirectProcess.cc.

236{
237 //Change the particle definition to the direct one
238 //------------------------------------------------
239 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
240 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
241
242 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
243
245 theDynPart->SetDefinition(theFwdParticleDef);
246
247
248 //Call the direct process
249 //----------------------
250 G4VParticleChange* partChange =theDirectProcess->AtRestDoIt( track, stepData );
251
252 //Restore the adjoint particle definition to the direct one
253 //------------------------------------------------
254 theDynPart->SetDefinition(adjPartDef);
255 theDynPart->SetPreAssignedDecayProducts(decayProducts);
256
257 return partChange;
258
259
260}
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0

◆ AtRestGetPhysicalInteractionLength()

G4double G4AdjointProcessEquivalentToDirectProcess::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 103 of file G4AdjointProcessEquivalentToDirectProcess.cc.

106{ //Change the particle definition to the direct one
107 //------------------------------------------------
108 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
109 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
110
111 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
113 theDynPart->SetDefinition(theFwdParticleDef);
114
115
116 //Call the direct process
117 //----------------------
118
119
120 G4double GPIL = theDirectProcess->AtRestGetPhysicalInteractionLength( track, condition );
121
122 //Restore the adjoint particle definition to the direct one
123 //------------------------------------------------
124 theDynPart->SetDefinition(adjPartDef);
125 theDynPart->SetPreAssignedDecayProducts(decayProducts);
126
127 return GPIL;
128
129
130}
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)=0

◆ BuildPhysicsTable()

void G4AdjointProcessEquivalentToDirectProcess::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 267 of file G4AdjointProcessEquivalentToDirectProcess.cc.

268{
269 return theDirectProcess->BuildPhysicsTable(*theFwdParticleDef);
270}
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:210

◆ EndTracking()

void G4AdjointProcessEquivalentToDirectProcess::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 316 of file G4AdjointProcessEquivalentToDirectProcess.cc.

317{
318 theDirectProcess->EndTracking();
319}
virtual void EndTracking()
Definition: G4VProcess.cc:137

◆ IsApplicable()

G4bool G4AdjointProcessEquivalentToDirectProcess::IsApplicable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 262 of file G4AdjointProcessEquivalentToDirectProcess.cc.

263{
264 return theDirectProcess->IsApplicable(*theFwdParticleDef);
265}
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:205

◆ PostStepDoIt()

G4VParticleChange * G4AdjointProcessEquivalentToDirectProcess::PostStepDoIt ( const G4Track track,
const G4Step stepData 
)
virtual

Implements G4VProcess.

Definition at line 177 of file G4AdjointProcessEquivalentToDirectProcess.cc.

179{
180 //Change the particle definition to the direct one
181 //------------------------------------------------
182 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
183 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
184
185 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
186
188 theDynPart->SetDefinition(theFwdParticleDef);
189
190
191 //Call the direct process
192 //----------------------
193
194 G4VParticleChange* partChange = theDirectProcess->PostStepDoIt( track, stepData );
195
196
197 //Restore the adjoint particle definition to the direct one
198 //------------------------------------------------
199 theDynPart->SetDefinition(adjPartDef);
200 theDynPart->SetPreAssignedDecayProducts(decayProducts);
201
202 return partChange;
203
204
205
206}
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0

◆ PostStepGetPhysicalInteractionLength()

G4double G4AdjointProcessEquivalentToDirectProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 132 of file G4AdjointProcessEquivalentToDirectProcess.cc.

136{
137 //Change the particle definition to the direct one
138 //------------------------------------------------
139 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
140 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
141
142 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
143
145 theDynPart->SetDefinition(theFwdParticleDef);
146
147
148 //Call the direct process
149 //----------------------
150
151
152 G4double GPIL = theDirectProcess->PostStepGetPhysicalInteractionLength( track,
153 previousStepSize,
154 condition );
155
156 //Restore the adjoint particle definition to the direct one
157 //------------------------------------------------
158 theDynPart->SetDefinition(adjPartDef);
159 theDynPart->SetPreAssignedDecayProducts(decayProducts);
160
161 return GPIL;
162
163
164}
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0

◆ PreparePhysicsTable()

void G4AdjointProcessEquivalentToDirectProcess::PreparePhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 272 of file G4AdjointProcessEquivalentToDirectProcess.cc.

273{
274 return theDirectProcess->PreparePhysicsTable(*theFwdParticleDef);
275}
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:217

◆ ResetNumberOfInteractionLengthLeft()

void G4AdjointProcessEquivalentToDirectProcess::ResetNumberOfInteractionLengthLeft ( )
virtual

Reimplemented from G4VProcess.

Definition at line 59 of file G4AdjointProcessEquivalentToDirectProcess.cc.

60{
61 theDirectProcess->ResetNumberOfInteractionLengthLeft();
62}
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:92

◆ RetrievePhysicsTable()

G4bool G4AdjointProcessEquivalentToDirectProcess::RetrievePhysicsTable ( const G4ParticleDefinition ,
const G4String directory,
G4bool  ascii = false 
)
virtual

Reimplemented from G4VProcess.

Definition at line 285 of file G4AdjointProcessEquivalentToDirectProcess.cc.

289{
290 return theDirectProcess->RetrievePhysicsTable(theFwdParticleDef, directory, ascii);
291}
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:236

◆ StartTracking()

void G4AdjointProcessEquivalentToDirectProcess::StartTracking ( G4Track track)
virtual

Reimplemented from G4VProcess.

Definition at line 293 of file G4AdjointProcessEquivalentToDirectProcess.cc.

294{
295 //Change the particle definition to the direct one
296 //------------------------------------------------
297 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track->GetDynamicParticle());
298 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
299
300 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
302 theDynPart->SetDefinition(theFwdParticleDef);
303
304 theDirectProcess->StartTracking(track);
305
306 //Restore the adjoint particle definition to the direct one
307 //------------------------------------------------
308 theDynPart->SetDefinition(adjPartDef);
309 theDynPart->SetPreAssignedDecayProducts(decayProducts);
310
311
312 return;
313
314}
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:125

◆ StorePhysicsTable()

G4bool G4AdjointProcessEquivalentToDirectProcess::StorePhysicsTable ( const G4ParticleDefinition ,
const G4String directory,
G4bool  ascii = false 
)
virtual

Reimplemented from G4VProcess.

Definition at line 277 of file G4AdjointProcessEquivalentToDirectProcess.cc.

281{
282 return theDirectProcess->StorePhysicsTable(theFwdParticleDef, directory, ascii);
283}
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:231

The documentation for this class was generated from the following files: