Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4WeightCutOffProcess Class Reference

#include <G4WeightCutOffProcess.hh>

Inheritance diagram for G4WeightCutOffProcess:
Collaboration diagram for G4WeightCutOffProcess:

Public Member Functions

 G4WeightCutOffProcess (G4double wsurvival, G4double wlimit, G4double isource, G4VIStore *istore, const G4String &aName="WeightCutOffProcess", G4bool para=false)
 
virtual ~G4WeightCutOffProcess ()
 
void SetParallelWorld (const G4String &parallelWorldName)
 
void SetParallelWorld (G4VPhysicalVolume *parallelWorld)
 
void StartTracking (G4Track *)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
const G4StringGetName () const
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- 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
 
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 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
 
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 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 58 of file G4WeightCutOffProcess.hh.

Constructor & Destructor Documentation

G4WeightCutOffProcess::G4WeightCutOffProcess ( G4double  wsurvival,
G4double  wlimit,
G4double  isource,
G4VIStore istore,
const G4String aName = "WeightCutOffProcess",
G4bool  para = false 
)

Definition at line 55 of file G4WeightCutOffProcess.cc.

61  : G4VProcess(aName),
62  fParticleChange(new G4ParticleChange),
63  fWeightSurvival(wsurvival),
64  fWeightLimit(wlimit),
65  fSourceImportance(isource),
66  fIStore(istore),
67  // fGCellFinder(aGCellFinder),
68  fGhostWorldName("NoParallelWorld"), fGhostWorld(0),
69  fGhostNavigator(0), fNavigatorID(-1), fFieldTrack('0'),
70  fParaflag(para), fEndTrack('0'), feLimited(kDoNot)
71 {
72  if (!fParticleChange)
73  {
74  G4Exception("G4WeightCutOffProcess::G4WeightCutOffProcess()",
75  "FatalError", FatalException,
76  "Failed to allocate G4ParticleChange !");
77  }
78 
79  G4VProcess::pParticleChange = fParticleChange;
80 
81  fGhostStep = new G4Step();
82  fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
83  fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
84 
85  fTransportationManager = G4TransportationManager::GetTransportationManager();
86  fPathFinder = G4PathFinder::GetInstance();
87 
88  if (verboseLevel>0)
89  {
90  G4cout << GetProcessName() << " is created " << G4endl;
91  }
92 }
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
G4StepPoint * GetPreStepPoint() const
G4GLOB_DLL std::ostream G4cout
Definition: G4Step.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
G4StepPoint * GetPostStepPoint() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4WeightCutOffProcess::~G4WeightCutOffProcess ( )
virtual

Definition at line 94 of file G4WeightCutOffProcess.cc.

95 {
96  delete fParticleChange;
97  // delete fGhostStep;
98 }

Member Function Documentation

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

Implements G4VProcess.

Definition at line 384 of file G4WeightCutOffProcess.cc.

385 {
386  // Dummy ParticleChange ie: does nothing
387  // Expecting G4Transportation to move the track
388  pParticleChange->Initialize(track);
389  return pParticleChange;
390 
391  // return 0;
392 }
virtual void Initialize(const G4Track &)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

Here is the call graph for this function:

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

Implements G4VProcess.

Definition at line 297 of file G4WeightCutOffProcess.cc.

300 {
301  if(fParaflag) {
302 
303  *selection = NotCandidateForSelection;
304  G4double returnedStep = DBL_MAX;
305 
306  if (previousStepSize > 0.)
307  { fGhostSafety -= previousStepSize; }
308  // else
309  // { fGhostSafety = -1.; }
310  if (fGhostSafety < 0.) fGhostSafety = 0.0;
311 
312  // ------------------------------------------
313  // Determination of the proposed STEP LENGTH:
314  // ------------------------------------------
315  if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
316  {
317  // I have no chance to limit
318  returnedStep = currentMinimumStep;
319  fOnBoundary = false;
320  proposedSafety = fGhostSafety - currentMinimumStep;
321  }
322  else // (currentMinimumStep > fGhostSafety: I may limit the Step)
323  {
324  G4FieldTrackUpdator::Update(&fFieldTrack,&track);
325  //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
326  // ComputeStep
327  //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
328  returnedStep
329  = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
330  track.GetCurrentStepNumber(),fGhostSafety,feLimited,
331  fEndTrack,track.GetVolume());
332  //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
333  if(feLimited == kDoNot)
334  {
335  // Track is not on the boundary
336  fOnBoundary = false;
337  fGhostSafety = fGhostNavigator->ComputeSafety(fEndTrack.GetPosition());
338  }
339  else
340  {
341  // Track is on the boundary
342  fOnBoundary = true;
343  proposedSafety = fGhostSafety;
344  }
345  //xbug? proposedSafety = fGhostSafety;
346  if(feLimited == kUnique || feLimited == kSharedOther) {
347  *selection = CandidateForSelection;
348  }else if (feLimited == kSharedTransport) {
349  returnedStep *= (1.0 + 1.0e-9);
350  // Expand to disable its selection in Step Manager comparison
351  }
352 
353  }
354 
355  // ----------------------------------------------
356  // Returns the fGhostSafety as the proposedSafety
357  // The SteppingManager will take care of keeping
358  // the smallest one.
359  // ----------------------------------------------
360  return returnedStep;
361 
362  } else {
363  return DBL_MAX;
364  //not sensible! return -1.0;
365  }
366 
367 }
G4int GetCurrentStepNumber() const
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4VPhysicalVolume * GetVolume() const
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
#define DBL_MAX
Definition: templates.hh:83
static void Update(G4FieldTrack *, const G4Track *)

Here is the call graph for this function:

G4VParticleChange * G4WeightCutOffProcess::AtRestDoIt ( const G4Track ,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 378 of file G4WeightCutOffProcess.cc.

379 {
380  return 0;
381 }
G4double G4WeightCutOffProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
virtual

Implements G4VProcess.

Definition at line 371 of file G4WeightCutOffProcess.cc.

373 {
374  return -1.0;
375 }
const G4String & G4WeightCutOffProcess::GetName ( ) const

Definition at line 291 of file G4WeightCutOffProcess.cc.

292 {
293  return theProcessName;
294 }
G4String theProcessName
Definition: G4VProcess.hh:335
G4VParticleChange * G4WeightCutOffProcess::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Implements G4VProcess.

Definition at line 188 of file G4WeightCutOffProcess.cc.

190 {
191  fParticleChange->Initialize(aTrack);
192 
193  if(fParaflag) {
194  fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
195  //xbug? fOnBoundary = false;
196  CopyStep(aStep);
197 
198  if(fOnBoundary)
199  {
200 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
201 // Locate the point and get new touchable
202 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
203  //?? fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
204  //?? step.GetPostStepPoint()->GetMomentumDirection());
205  fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
206 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
207  }
208  else
209  {
210  // Do I need this ??????????????????????????????????????????????????????????
211  // fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
212  // ?????????????????????????????????????????????????????????????????????????
213 
214  // fPathFinder->ReLocate(track.GetPosition());
215 
216  // reuse the touchable
217  fNewGhostTouchable = fOldGhostTouchable;
218  }
219 
220  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
221  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
222 
223  }
224 
225  if(fParaflag) {
226  G4GeometryCell postCell(*(fGhostPostStepPoint->GetPhysicalVolume()),
227  fGhostPostStepPoint->GetTouchable()->GetReplicaNumber());
228 
229 
230  // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(aStep);
231  // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(fGhostStep);
232  G4double R = fSourceImportance;
233  if (fIStore)
234  {
235  G4double i = fIStore->GetImportance(postCell);
236  if (i>0)
237  {
238  R/=i;
239  }
240  }
241  G4double w = aTrack.GetWeight();
242  if (w<R*fWeightLimit)
243  {
244  G4double ws = fWeightSurvival*R;
245  G4double p = w/(ws);
246  if (G4UniformRand()<p)
247  {
248  fParticleChange->ProposeTrackStatus(fStopAndKill);
249  }
250  else
251  {
252  fParticleChange->ProposeWeight(ws);
253  }
254  }
255  } else {
256 
257  G4GeometryCell postCell(*(aStep.GetPostStepPoint()->GetPhysicalVolume()),
259 
260  // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(aStep);
261  // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(fGhostStep);
262  G4double R = fSourceImportance;
263  if (fIStore)
264  {
265  G4double i = fIStore->GetImportance(postCell);
266  if (i>0)
267  {
268  R/=i;
269  }
270  }
271  G4double w = aTrack.GetWeight();
272  if (w<R*fWeightLimit)
273  {
274  G4double ws = fWeightSurvival*R;
275  G4double p = w/(ws);
276  if (G4UniformRand()<p)
277  {
278  fParticleChange->ProposeTrackStatus(fStopAndKill);
279  }
280  else
281  {
282  fParticleChange->ProposeWeight(ws);
283  }
284  }
285  }
286 
287  return fParticleChange;
288 
289 }
const char * p
Definition: xmltok.h:285
G4TouchableHandle CreateTouchableHandle(G4int navId) const
const G4VTouchable * GetTouchable() const
void ProposeWeight(G4double finalWeight)
#define G4UniformRand()
Definition: Randomize.hh:97
G4VPhysicalVolume * GetPhysicalVolume() const
virtual G4double GetImportance(const G4GeometryCell &gCell) const =0
virtual void Initialize(const G4Track &)
G4StepPoint * GetPostStepPoint() const
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58
G4double GetWeight() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const

Here is the call graph for this function:

G4double G4WeightCutOffProcess::PostStepGetPhysicalInteractionLength ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 176 of file G4WeightCutOffProcess.cc.

178 {
179 // *condition = Forced;
180 // return kInfinity;
181 
182 // *condition = StronglyForced;
183  *condition = Forced;
184  return DBL_MAX;
185 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
void G4WeightCutOffProcess::SetParallelWorld ( const G4String parallelWorldName)

Definition at line 107 of file G4WeightCutOffProcess.cc.

108 {
109 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
110 // Get pointers of the parallel world and its navigator
111 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
112  fGhostWorldName = parallelWorldName;
113  fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
114  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
115 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
116 }
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
G4Navigator * GetNavigator(const G4String &worldName)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4WeightCutOffProcess::SetParallelWorld ( G4VPhysicalVolume parallelWorld)

Definition at line 119 of file G4WeightCutOffProcess.cc.

120 {
121 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
122 // Get pointer of navigator
123 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
124  fGhostWorldName = parallelWorld->GetName();
125  fGhostWorld = parallelWorld;
126  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
127 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
128 }
const G4String & GetName() const
G4Navigator * GetNavigator(const G4String &worldName)

Here is the call graph for this function:

void G4WeightCutOffProcess::StartTracking ( G4Track trk)
virtual

Reimplemented from G4VProcess.

Definition at line 135 of file G4WeightCutOffProcess.cc.

136 {
137 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
138 // Activate navigator and get the navigator ID
139 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
140 // G4cout << " G4ParallelWorldScoringProcess::StartTracking" << G4endl;
141 
142  if(fParaflag) {
143  if(fGhostNavigator)
144  { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
145  else
146  {
147  G4Exception("G4WeightCutOffProcess::StartTracking",
148  "ProcParaWorld000",FatalException,
149  "G4WeightCutOffProcess is used for tracking without having a parallel world assigned");
150  }
151 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
152 
153 // G4cout << "G4ParallelWorldScoringProcess::StartTracking <<<<<<<<<<<<<<<<<< " << G4endl;
154 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
155 // Let PathFinder initialize
156 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
157  fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
158 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
159 
160 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
161 // Setup initial touchables for the first step
162 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
163  fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
164  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
165  fNewGhostTouchable = fOldGhostTouchable;
166  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
167 
168  // Initialize
169  fGhostSafety = -1.;
170  fOnBoundary = false;
171  }
172 }
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=0)
const G4ThreeVector & GetPosition() const
G4TouchableHandle CreateTouchableHandle(G4int navId) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int ActivateNavigator(G4Navigator *aNavigator)
const G4ThreeVector & GetMomentumDirection() const
void SetTouchableHandle(const G4TouchableHandle &apValue)

Here is the call graph for this function:


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