Geant4  10.03.p01
 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 // $Id: G4AdjointProcessEquivalentToDirectProcess.cc 66892 2013-01-17 10:57:59Z gunter $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 class implementation file
32 //
33 // Class Description
34 //
35 // This class is for adjoint process equivalent to direct process
36 
37 // ------------------------------------------------------------
38 // Created by L.Desorgher 25 Sept. 2009 Inspired from G4WrapperProcess
39 // ------------------------------------------------------------
40 
42 #include "G4DynamicParticle.hh"
44  G4VProcess* aProcess,
45  G4ParticleDefinition* fwd_particle_def)
46 :G4VProcess(aName)
47 {
48  theDirectProcess =aProcess;
49  theProcessType = theDirectProcess->GetProcessType();
50  theFwdParticleDef = fwd_particle_def;
51 }
52 
53 
55 {
56  if (theDirectProcess!=0) delete theDirectProcess;
57 }
58 
60 {
61  theDirectProcess->ResetNumberOfInteractionLengthLeft();
62 }
63 
66  G4double previousStepSize,
67  G4double currentMinimumStep,
68  G4double& proposedSafety,
69  G4GPILSelection* selection )
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 }
102 
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());
112  theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
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 }
131 
134  G4double previousStepSize,
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 
144  theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
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 }
165 /*
166 
167 void G4AdjointProcessEquivalentToDirectProcess::SetProcessManager(const G4ProcessManager* procMan)
168 {
169  theDirectProcess->SetProcessManager(procMan);
170 }
171 
172 const G4ProcessManager* G4AdjointProcessEquivalentToDirectProcess::GetProcessManager()
173 {
174  return theDirectProcess->GetProcessManager();
175 }
176 */
178  const G4Step& stepData )
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 
187  theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
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 }
207 
209  const G4Step& stepData )
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 
218  theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
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 }
233 
235  const G4Step& stepData )
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 
244  theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
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 }
261 
263 {
264  return theDirectProcess->IsApplicable(*theFwdParticleDef);
265 }
266 
268 {
269  return theDirectProcess->BuildPhysicsTable(*theFwdParticleDef);
270 }
271 
273 {
274  return theDirectProcess->PreparePhysicsTable(*theFwdParticleDef);
275 }
276 
279  const G4String& directory,
280  G4bool ascii)
281 {
282  return theDirectProcess->StorePhysicsTable(theFwdParticleDef, directory, ascii);
283 }
284 
287  const G4String& directory,
288  G4bool ascii)
289 {
290  return theDirectProcess->RetrievePhysicsTable(theFwdParticleDef, directory, ascii);
291 }
292 
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());
301  theDynPart->SetPreAssignedDecayProducts((G4DecayProducts*)(0));
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 }
315 
317 {
318  theDirectProcess->EndTracking();
319 }
320 
G4double condition(const G4ErrorSymMatrix &m)
G4ProcessType theProcessType
Definition: G4VProcess.hh:340
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)=0
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
const G4DynamicParticle * GetDynamicParticle() const
const G4DecayProducts * GetPreAssignedDecayProducts() const
G4ParticleDefinition * GetDefinition() const
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:101
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:217
bool G4bool
Definition: G4Types.hh:79
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)
Definition: G4Step.hh:76
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:236
G4AdjointProcessEquivalentToDirectProcess(const G4String &aName, G4VProcess *aProcess, G4ParticleDefinition *fwd_particle_def)
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:210
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0
virtual G4bool IsApplicable(const G4ParticleDefinition &)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:231
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
virtual void EndTracking()
Definition: G4VProcess.cc:113
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
void SetPreAssignedDecayProducts(G4DecayProducts *aDecayProducts)
G4ForceCondition
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:205
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4GPILSelection
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0