Geant4 10.7.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 ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool 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
 
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 54 of file G4AdjointProcessEquivalentToDirectProcess.hh.

Constructor & Destructor Documentation

◆ G4AdjointProcessEquivalentToDirectProcess()

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

Definition at line 42 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

◆ ~G4AdjointProcessEquivalentToDirectProcess()

G4AdjointProcessEquivalentToDirectProcess::~G4AdjointProcessEquivalentToDirectProcess ( )
virtual

Definition at line 53 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 207 of file G4AdjointProcessEquivalentToDirectProcess.cc.

209{
210 //Change the particle definition to the direct one
211 //------------------------------------------------
212 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
213 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
214
215 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
216
218 theDynPart->SetDefinition(theFwdParticleDef);
219
220
221 //Call the direct process
222 //----------------------
223 G4VParticleChange* partChange =theDirectProcess->AlongStepDoIt( track, stepData );
224
225 //Restore the adjoint particle definition to the direct one
226 //------------------------------------------------
227 theDynPart->SetDefinition(adjPartDef);
228 theDynPart->SetPreAssignedDecayProducts(decayProducts);
229
230 return partChange;
231}
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 63 of file G4AdjointProcessEquivalentToDirectProcess.cc.

69{
70
71
72 //Change the particle definition to the direct one
73 //------------------------------------------------
74 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
75 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
76
77 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
79 theDynPart->SetDefinition(theFwdParticleDef);
80
81
82 //Call the direct process
83 //----------------------
84 G4double GPIL = theDirectProcess->
86 previousStepSize,
87 currentMinimumStep,
88 proposedSafety,
89 selection );
90
91
92 //Restore the adjoint particle definition to the direct one
93 //------------------------------------------------
94 theDynPart->SetDefinition(adjPartDef);
95 theDynPart->SetPreAssignedDecayProducts(decayProducts);
96
97
98 return GPIL;
99
100}
double G4double
Definition: G4Types.hh:83
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 233 of file G4AdjointProcessEquivalentToDirectProcess.cc.

235{
236 //Change the particle definition to the direct one
237 //------------------------------------------------
238 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
239 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
240
241 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
242
244 theDynPart->SetDefinition(theFwdParticleDef);
245
246
247 //Call the direct process
248 //----------------------
249 G4VParticleChange* partChange =theDirectProcess->AtRestDoIt( track, stepData );
250
251 //Restore the adjoint particle definition to the direct one
252 //------------------------------------------------
253 theDynPart->SetDefinition(adjPartDef);
254 theDynPart->SetPreAssignedDecayProducts(decayProducts);
255
256 return partChange;
257
258
259}
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 102 of file G4AdjointProcessEquivalentToDirectProcess.cc.

105{ //Change the particle definition to the direct one
106 //------------------------------------------------
107 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
108 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
109
110 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
112 theDynPart->SetDefinition(theFwdParticleDef);
113
114
115 //Call the direct process
116 //----------------------
117
118
119 G4double GPIL = theDirectProcess->AtRestGetPhysicalInteractionLength( track, condition );
120
121 //Restore the adjoint particle definition to the direct one
122 //------------------------------------------------
123 theDynPart->SetDefinition(adjPartDef);
124 theDynPart->SetPreAssignedDecayProducts(decayProducts);
125
126 return GPIL;
127
128
129}
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 266 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

◆ EndTracking()

void G4AdjointProcessEquivalentToDirectProcess::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 315 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

◆ IsApplicable()

G4bool G4AdjointProcessEquivalentToDirectProcess::IsApplicable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 261 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 176 of file G4AdjointProcessEquivalentToDirectProcess.cc.

178{
179 //Change the particle definition to the direct one
180 //------------------------------------------------
181 G4DynamicParticle* theDynPart = const_cast<G4DynamicParticle*> (track.GetDynamicParticle());
182 G4ParticleDefinition* adjPartDef = theDynPart->GetDefinition();
183
184 G4DecayProducts* decayProducts = const_cast<G4DecayProducts*> (theDynPart->GetPreAssignedDecayProducts());
185
187 theDynPart->SetDefinition(theFwdParticleDef);
188
189
190 //Call the direct process
191 //----------------------
192
193 G4VParticleChange* partChange = theDirectProcess->PostStepDoIt( track, stepData );
194
195
196 //Restore the adjoint particle definition to the direct one
197 //------------------------------------------------
198 theDynPart->SetDefinition(adjPartDef);
199 theDynPart->SetPreAssignedDecayProducts(decayProducts);
200
201 return partChange;
202
203
204
205}
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 131 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

◆ PreparePhysicsTable()

void G4AdjointProcessEquivalentToDirectProcess::PreparePhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 271 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

◆ ResetNumberOfInteractionLengthLeft()

void G4AdjointProcessEquivalentToDirectProcess::ResetNumberOfInteractionLengthLeft ( )
virtual

Reimplemented from G4VProcess.

Definition at line 58 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

◆ RetrievePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 284 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

◆ StartTracking()

void G4AdjointProcessEquivalentToDirectProcess::StartTracking ( G4Track track)
virtual

Reimplemented from G4VProcess.

Definition at line 292 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

◆ StorePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 276 of file G4AdjointProcessEquivalentToDirectProcess.cc.

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

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