Geant4 11.2.2
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)
 
 ~G4AdjointProcessEquivalentToDirectProcess () override
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData) override
 
G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
void StartTracking (G4Track *) override
 
void EndTracking () override
 
void ResetNumberOfInteractionLengthLeft () override
 
 G4AdjointProcessEquivalentToDirectProcess (G4AdjointProcessEquivalentToDirectProcess &)=delete
 
G4AdjointProcessEquivalentToDirectProcessoperator= (const G4AdjointProcessEquivalentToDirectProcess &right)=delete
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
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)
 
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 const G4VProcessGetCreatorProcess () const
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

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

Detailed Description

Definition at line 41 of file G4AdjointProcessEquivalentToDirectProcess.hh.

Constructor & Destructor Documentation

◆ G4AdjointProcessEquivalentToDirectProcess() [1/2]

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

Definition at line 40 of file G4AdjointProcessEquivalentToDirectProcess.cc.

44 : G4VProcess(aName)
45{
46 fDirectProcess = aProcess;
47 theProcessType = fDirectProcess->GetProcessType();
48 fFwdParticleDef = fwd_particle_def;
49}
G4ProcessType theProcessType
G4ProcessType GetProcessType() const

◆ ~G4AdjointProcessEquivalentToDirectProcess()

G4AdjointProcessEquivalentToDirectProcess::~G4AdjointProcessEquivalentToDirectProcess ( )
override

Definition at line 51 of file G4AdjointProcessEquivalentToDirectProcess.cc.

53{
54 if(fDirectProcess != nullptr)
55 delete fDirectProcess;
56}

◆ G4AdjointProcessEquivalentToDirectProcess() [2/2]

G4AdjointProcessEquivalentToDirectProcess::G4AdjointProcessEquivalentToDirectProcess ( G4AdjointProcessEquivalentToDirectProcess & )
delete

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 167 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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}
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 )
overridevirtual

Implements G4VProcess.

Definition at line 64 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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}
double G4double
Definition G4Types.hh:83
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0

◆ AtRestDoIt()

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

Implements G4VProcess.

Definition at line 192 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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}
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0

◆ AtRestGetPhysicalInteractionLength()

G4double G4AdjointProcessEquivalentToDirectProcess::AtRestGetPhysicalInteractionLength ( const G4Track & track,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 93 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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}
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)=0

◆ BuildPhysicsTable()

void G4AdjointProcessEquivalentToDirectProcess::BuildPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 222 of file G4AdjointProcessEquivalentToDirectProcess.cc.

224{
225 return fDirectProcess->BuildPhysicsTable(*fFwdParticleDef);
226}
virtual void BuildPhysicsTable(const G4ParticleDefinition &)

◆ EndTracking()

void G4AdjointProcessEquivalentToDirectProcess::EndTracking ( )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 268 of file G4AdjointProcessEquivalentToDirectProcess.cc.

269{
270 fDirectProcess->EndTracking();
271}
virtual void EndTracking()

◆ IsApplicable()

G4bool G4AdjointProcessEquivalentToDirectProcess::IsApplicable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 216 of file G4AdjointProcessEquivalentToDirectProcess.cc.

218{
219 return fDirectProcess->IsApplicable(*fFwdParticleDef);
220}
virtual G4bool IsApplicable(const G4ParticleDefinition &)

◆ operator=()

G4AdjointProcessEquivalentToDirectProcess & G4AdjointProcessEquivalentToDirectProcess::operator= ( const G4AdjointProcessEquivalentToDirectProcess & right)
delete

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 143 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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}
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0

◆ PostStepGetPhysicalInteractionLength()

G4double G4AdjointProcessEquivalentToDirectProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 118 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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}
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0

◆ PreparePhysicsTable()

void G4AdjointProcessEquivalentToDirectProcess::PreparePhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 228 of file G4AdjointProcessEquivalentToDirectProcess.cc.

230{
231 return fDirectProcess->PreparePhysicsTable(*fFwdParticleDef);
232}
virtual void PreparePhysicsTable(const G4ParticleDefinition &)

◆ ResetNumberOfInteractionLengthLeft()

void G4AdjointProcessEquivalentToDirectProcess::ResetNumberOfInteractionLengthLeft ( )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 58 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

◆ RetrievePhysicsTable()

G4bool G4AdjointProcessEquivalentToDirectProcess::RetrievePhysicsTable ( const G4ParticleDefinition * ,
const G4String & directory,
G4bool ascii = false )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 240 of file G4AdjointProcessEquivalentToDirectProcess.cc.

242{
243 return fDirectProcess->RetrievePhysicsTable(fFwdParticleDef, directory,
244 ascii);
245}
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)

◆ StartTracking()

void G4AdjointProcessEquivalentToDirectProcess::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 247 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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}
virtual void StartTracking(G4Track *)
Definition G4VProcess.cc:87

◆ StorePhysicsTable()

G4bool G4AdjointProcessEquivalentToDirectProcess::StorePhysicsTable ( const G4ParticleDefinition * ,
const G4String & directory,
G4bool ascii = false )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 234 of file G4AdjointProcessEquivalentToDirectProcess.cc.

236{
237 return fDirectProcess->StorePhysicsTable(fFwdParticleDef, directory, ascii);
238}
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)

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