source: trunk/source/processes/biasing/src/G4WeightCutOffProcess.cc @ 846

Last change on this file since 846 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 12.8 KB
Line 
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: G4WeightCutOffProcess.cc,v 1.2 2007/06/01 09:16:34 ahoward Exp $
28// GEANT4 tag $Name:  $
29//
30// ----------------------------------------------------------------------
31// GEANT 4 class source file
32//
33// G4WeightCutOffProcess.cc
34//
35// ----------------------------------------------------------------------
36
37#include "G4WeightCutOffProcess.hh"
38//#include "G4VScorer.hh"
39#include "G4GeometryCellStep.hh"
40#include "G4GCellFinder.hh"
41#include "G4TouchableHandle.hh"
42#include "G4VIStore.hh"
43
44#include "G4Step.hh"
45#include "G4Navigator.hh"
46#include "G4VTouchable.hh"
47#include "G4VPhysicalVolume.hh"
48#include "G4ParticleChange.hh"
49#include "G4PathFinder.hh"
50#include "G4TransportationManager.hh"
51#include "G4StepPoint.hh"
52#include "G4FieldTrackUpdator.hh"
53
54
55G4WeightCutOffProcess::
56G4WeightCutOffProcess(G4double wsurvival,
57                      G4double wlimit,
58                      G4double isource,
59                      G4VIStore *istore,
60                      const G4VGCellFinder &aGCellFinder,
61                      const G4String &aName, G4bool para)
62  : G4VProcess(aName), 
63    fParticleChange(new G4ParticleChange),
64    fWeightSurvival(wsurvival),
65    fWeightLimit(wlimit),
66    fSourceImportance(isource),
67    fIStore(istore),
68    fGCellFinder(aGCellFinder),
69   fGhostNavigator(0), fNavigatorID(-1), fFieldTrack('0')
70{
71  if (!fParticleChange)
72  {
73    G4Exception("G4WeightCutOffProcess::G4WeightCutOffProcess()",
74                "FatalError", FatalException,
75                "Failed to allocate G4ParticleChange !");
76  }
77
78  G4VProcess::pParticleChange = fParticleChange;
79
80  fGhostStep = new G4Step();
81  fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
82  fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
83
84  fTransportationManager = G4TransportationManager::GetTransportationManager();
85  fPathFinder = G4PathFinder::GetInstance();
86
87  if (verboseLevel>0)
88  {
89    G4cout << GetProcessName() << " is created " << G4endl;
90  }
91
92  paraflag = para;
93
94}
95
96G4WeightCutOffProcess::~G4WeightCutOffProcess()
97{
98  delete fParticleChange;
99  delete fGhostStep;
100}
101
102
103//------------------------------------------------------
104//
105// SetParallelWorld
106//
107//------------------------------------------------------
108void G4WeightCutOffProcess::
109SetParallelWorld(G4String parallelWorldName)
110{
111//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
112// Get pointers of the parallel world and its navigator
113//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
114  fGhostWorldName = parallelWorldName;
115  fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
116  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
117//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
118}
119
120void G4WeightCutOffProcess::
121SetParallelWorld(G4VPhysicalVolume* parallelWorld)
122{
123//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
124// Get pointer of navigator
125//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
126  fGhostWorldName = parallelWorld->GetName();
127  fGhostWorld = parallelWorld;
128  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
129//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
130}
131
132//------------------------------------------------------
133//
134// StartTracking
135//
136//------------------------------------------------------
137void G4WeightCutOffProcess::StartTracking(G4Track* trk)
138{
139//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
140// Activate navigator and get the navigator ID
141//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
142// G4cout << " G4ParallelWorldScoringProcess::StartTracking" << G4endl;
143
144  if(paraflag) {
145    if(fGhostNavigator)
146      { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
147    else
148      {
149        G4Exception("G4WeightCutOffProcess::StartTracking",
150                    "ProcParaWorld000",FatalException,
151                    "G4WeightCutOffProcess is used for tracking without having a parallel world assigned");
152      }
153//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
154
155// G4cout << "G4ParallelWorldScoringProcess::StartTracking <<<<<<<<<<<<<<<<<< " << G4endl;
156//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
157// Let PathFinder initialize
158//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
159    fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
160//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
161
162//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
163// Setup initial touchables for the first step
164//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
165    fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
166    fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
167    fNewGhostTouchable = fOldGhostTouchable;
168    fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
169   
170    // Initialize
171    fGhostSafety = -1.;
172    fOnBoundary = false;
173  }
174}
175
176
177G4double G4WeightCutOffProcess::
178PostStepGetPhysicalInteractionLength(const G4Track &,
179                                     G4double, G4ForceCondition* condition)
180{
181//   *condition = Forced;
182//   return kInfinity;
183
184//  *condition = StronglyForced;
185  *condition = Forced;
186  return DBL_MAX;
187}
188 
189G4VParticleChange * 
190G4WeightCutOffProcess::PostStepDoIt(const G4Track& aTrack, 
191                                    const G4Step &aStep)
192{
193  fParticleChange->Initialize(aTrack);
194
195  if(paraflag) {
196    fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
197    //xbug?    fOnBoundary = false;
198    CopyStep(aStep);
199
200    if(fOnBoundary)
201      {
202//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
203// Locate the point and get new touchable
204//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
205  //??  fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
206  //??                      step.GetPostStepPoint()->GetMomentumDirection());
207        fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
208//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
209      }
210    else
211      {
212        // Do I need this ??????????????????????????????????????????????????????????
213        // fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
214        // ?????????????????????????????????????????????????????????????????????????
215       
216        // fPathFinder->ReLocate(track.GetPosition());
217       
218        // reuse the touchable
219        fNewGhostTouchable = fOldGhostTouchable;
220      }
221   
222    fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
223    fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
224
225  }
226
227  G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(aStep);
228  //  G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(fGhostStep);
229  G4double R = fSourceImportance;
230  if (fIStore)
231  {
232    G4double i = fIStore->GetImportance(postCell);
233    if (i>0)
234    {
235      R/=i;
236    }
237  }
238  G4double w = aTrack.GetWeight();
239  if (w<R*fWeightLimit)
240  {
241    G4double ws = fWeightSurvival*R;
242    G4double p = w/(ws);
243    if (G4UniformRand()<p)
244    {
245      fParticleChange->ProposeTrackStatus(fStopAndKill);
246    }
247    else
248    {
249      fParticleChange->ProposeWeight(ws);
250    }                 
251  }
252  return fParticleChange;
253}
254
255const G4String &G4WeightCutOffProcess::GetName() const
256{
257  return theProcessName;
258}
259
260G4double G4WeightCutOffProcess::
261AlongStepGetPhysicalInteractionLength(
262            const G4Track& track, G4double  previousStepSize, G4double  currentMinimumStep,
263            G4double& proposedSafety, G4GPILSelection* selection)
264{
265  if(paraflag) {
266    static G4FieldTrack endTrack('0');
267    static ELimited eLimited;
268   
269    *selection = NotCandidateForSelection;
270    G4double returnedStep = DBL_MAX;
271   
272    if (previousStepSize > 0.)
273      { fGhostSafety -= previousStepSize; }
274    //  else
275    //  { fGhostSafety = -1.; }
276    if (fGhostSafety < 0.) fGhostSafety = 0.0;
277   
278    // ------------------------------------------
279    // Determination of the proposed STEP LENGTH:
280    // ------------------------------------------
281    if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
282      {
283        // I have no chance to limit
284        returnedStep = currentMinimumStep;
285        fOnBoundary = false;
286        proposedSafety = fGhostSafety - currentMinimumStep;
287      }
288    else // (currentMinimumStep > fGhostSafety: I may limit the Step)
289      {
290        G4FieldTrackUpdator::Update(&fFieldTrack,&track);
291        //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
292        // ComputeStep
293        //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
294        returnedStep
295          = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
296                                     track.GetCurrentStepNumber(),fGhostSafety,eLimited,
297                                     endTrack,track.GetVolume());
298        //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
299        if(eLimited == kDoNot)
300          {
301            // Track is not on the boundary
302            fOnBoundary = false;
303            fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
304          }
305        else
306          {
307            // Track is on the boundary
308            fOnBoundary = true;
309            proposedSafety = fGhostSafety;
310          }
311        //xbug? proposedSafety = fGhostSafety;
312        if(eLimited == kUnique || eLimited == kSharedOther) {
313          *selection = CandidateForSelection;
314        }else if (eLimited == kSharedTransport) { 
315          returnedStep *= (1.0 + 1.0e-9); 
316          // Expand to disable its selection in Step Manager comparison
317        }
318       
319      }
320
321  // ----------------------------------------------
322  // Returns the fGhostSafety as the proposedSafety
323  // The SteppingManager will take care of keeping
324  // the smallest one.
325  // ----------------------------------------------
326    return returnedStep;
327
328  } else {
329    return DBL_MAX;
330    //not sensible!    return -1.0;
331  }
332
333}
334
335 
336G4double G4WeightCutOffProcess::
337AtRestGetPhysicalInteractionLength(const G4Track& ,
338                                   G4ForceCondition*)
339{
340  return -1.0;
341}
342 
343G4VParticleChange*
344G4WeightCutOffProcess::AtRestDoIt(const G4Track&, const G4Step&)
345{
346  return 0;
347}
348
349G4VParticleChange*
350G4WeightCutOffProcess::AlongStepDoIt(const G4Track& track, const G4Step&)
351{
352  // Dummy ParticleChange ie: does nothing
353  // Expecting G4Transportation to move the track
354  pParticleChange->Initialize(track);
355  return pParticleChange; 
356
357  //  return 0;
358}
359
360void G4WeightCutOffProcess::CopyStep(const G4Step & step)
361{
362  fGhostStep->SetTrack(step.GetTrack());
363  fGhostStep->SetStepLength(step.GetStepLength());
364  fGhostStep->SetTotalEnergyDeposit(step.GetTotalEnergyDeposit());
365  fGhostStep->SetControlFlag(step.GetControlFlag());
366
367  *fGhostPreStepPoint = *(step.GetPreStepPoint());
368  *fGhostPostStepPoint = *(step.GetPostStepPoint());
369
370//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
371// Set StepStatus for ghost world
372//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
373  if(fOnBoundary)
374  { fGhostPostStepPoint->SetStepStatus(fGeomBoundary); }
375  else if(fGhostPostStepPoint->GetStepStatus()==fGeomBoundary)
376  { fGhostPostStepPoint->SetStepStatus(fPostStepDoItProc); }
377//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
378}
Note: See TracBrowser for help on using the repository browser.