Geant4  10.02
G4BOptrForceCollision.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 #include "G4BOptrForceCollision.hh"
29 #include "G4PhysicsModelCatalog.hh"
30 
34 #include "G4BOptnCloning.hh"
35 
36 #include "G4Step.hh"
37 #include "G4StepPoint.hh"
38 #include "G4VProcess.hh"
39 
40 #include "G4ParticleDefinition.hh"
41 #include "G4ParticleTable.hh"
42 
43 #include "G4SystemOfUnits.hh"
44 
45 // -- §§ consider calling other constructor, thanks to C++11
47  : G4VBiasingOperator(name),
48  fForceCollisionModelID(-1),
49  fCurrentTrackData(nullptr),
50  fSetup(true)
51 {
53  fCloningOperation = new G4BOptnCloning("Cloning");
55 
56  if ( fParticleToBias == 0 )
57  {
59  ed << " Particle `" << particleName << "' not found !" << G4endl;
60  G4Exception(" G4BOptrForceCollision::G4BOptrForceCollision(...)",
61  "BIAS.GEN.07",
63  ed);
64  }
65 }
66 
67 
69  : G4VBiasingOperator(name),
70  fForceCollisionModelID(-1),
71  fCurrentTrackData(nullptr),
72  fSetup(true)
73 {
75  fCloningOperation = new G4BOptnCloning("Cloning");
76  fParticleToBias = particle;
77 }
78 
79 
81 {
82  for ( std::map< const G4BiasingProcessInterface*, G4BOptnForceFreeFlight* >::iterator it = fFreeFlightOperations.begin() ;
83  it != fFreeFlightOperations.end() ;
84  it++ ) delete (*it).second;
86  delete fCloningOperation;
87 }
88 
89 
91 {
92  // -- Create ID for force collision:
94 
95  // -- build free flight operations:
97 }
98 
99 
101 {
102  // -- start by remembering processes under biasing, and create needed biasing operations:
103  if ( fSetup )
104  {
105  // -- get back ID created in master thread in case of MT (or reget just created ID in sequential mode):
106  fForceCollisionModelID = G4PhysicsModelCatalog::Register("GenBiasForceCollision");
107 
108  const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
109  const G4BiasingProcessSharedData* interfaceProcessSharedData = G4BiasingProcessInterface::GetSharedData( processManager );
110  if ( interfaceProcessSharedData ) // -- sharedData tested, as is can happen a user attaches an operator
111  // -- to a volume but without defining BiasingProcessInterface processes.
112  {
113  for ( size_t i = 0 ; i < (interfaceProcessSharedData->GetPhysicsBiasingProcessInterfaces()).size(); i++ )
114  {
115  const G4BiasingProcessInterface* wrapperProcess =
116  (interfaceProcessSharedData->GetPhysicsBiasingProcessInterfaces())[i];
117  G4String operationName = "FreeFlight-"+wrapperProcess->GetWrappedProcess()->GetProcessName();
118  fFreeFlightOperations[wrapperProcess] = new G4BOptnForceFreeFlight(operationName);
119  }
120  }
121  fSetup = false;
122  }
123 }
124 
125 
127 {
128 }
129 
130 
132 {
133  // -- does nothing if particle is not of requested type:
134  if ( track->GetDefinition() != fParticleToBias ) return 0;
135 
136  // -- trying to get auxiliary track data...
137  if ( fCurrentTrackData == nullptr )
138  {
139  // ... and if the track has no aux. track data, it means its biasing is not started yet (note that cloning is to happen first):
141  if ( fCurrentTrackData == nullptr ) return nullptr;
142  }
143 
144 
145  // -- Send force free flight to the callingProcess:
146  // ------------------------------------------------
147  // -- The track has been cloned in the previous step, it has now to be
148  // -- forced for a free flight.
149  // -- This track will fly with 0.0 weight during its forced flight:
150  // -- it is to forbid the double counting with the force interaction track.
151  // -- Its weight is restored at the end of its free flight, this weight
152  // -- being its initial weight * the weight for the free flight travel,
153  // -- this last one being per process. The initial weight is common, and is
154  // -- arbitrary asked to the first operation to take care of it.
156  {
157  G4BOptnForceFreeFlight* operation = fFreeFlightOperations[callingProcess];
158  if ( callingProcess->GetWrappedProcess()->GetCurrentInteractionLength() < DBL_MAX/10. )
159  {
160  // -- the initial track weight will be restored only by the first DoIt free flight:
162  return operation;
163  }
164  else
165  {
166  return nullptr;
167  }
168  }
169 
170 
171  // -- Send force interaction operation to the callingProcess:
172  // ----------------------------------------------------------
173  // -- at this level, a copy of the track entering the volume was
174  // -- generated (borned) earlier. This copy will make the forced
175  // -- interaction in the volume.
177  {
178  // -- Remember if this calling process is the first of the physics wrapper in
179  // -- the PostStepGPIL loop (using default argument of method below):
180  G4bool isFirstPhysGPIL = callingProcess-> GetIsFirstPostStepGPILInterface();
181 
182  // -- [*first process*] Initialize or update force interaction operation:
183  if ( isFirstPhysGPIL )
184  {
185  // -- first step of cloned track, initialize the forced interaction operation:
187  else
188  {
190  {
191  // -- means that some other physics process, not under control of the forced interaction operation,
192  // -- has occured, need to re-initialize the operation as distance to boundary has changed.
193  // -- [ Note the re-initialization is only possible for a Markovian law. ]
195  }
196  else
197  {
198  // -- means that some other non-physics process (biasing or not, like step limit), has occured,
199  // -- but track conserves its momentum direction, only need to reduced the maximum distance for
200  // -- forced interaction.
201  // -- [ Note the update is only possible for a Markovian law. ]
203  }
204  }
205  }
206 
207  // -- [*all processes*] Sanity check : it may happen in limit cases that distance to
208  // -- out is zero, weight would be infinite in this case: abort forced interaction
209  // -- and abandon biasing.
211  {
213  return 0;
214  }
215 
216  // -- [* first process*] collect cross-sections and sample force law to determine interaction length
217  // -- and winning process:
218  if ( isFirstPhysGPIL )
219  {
220  // -- collect cross-sections:
221  // -- ( Remember that the first of the G4BiasingProcessInterface triggers the update
222  // -- of these cross-sections )
223  const G4BiasingProcessSharedData* sharedData = callingProcess->GetSharedData();
224  for ( size_t i = 0 ; i < (sharedData->GetPhysicsBiasingProcessInterfaces()).size() ; i++ )
225  {
226  const G4BiasingProcessInterface* wrapperProcess = ( sharedData->GetPhysicsBiasingProcessInterfaces() )[i];
227  G4double interactionLength = wrapperProcess->GetWrappedProcess()->GetCurrentInteractionLength();
228  // -- keep only well defined cross-sections, other processes are ignored. These are not pathological cases
229  // -- but cases where a threhold effect par example (eg pair creation) may be at play. (**!**)
230  if ( interactionLength < DBL_MAX/10. )
231  fSharedForceInteractionOperation->AddCrossSection( wrapperProcess->GetWrappedProcess(), 1.0/interactionLength );
232  }
233  // -- sample the shared law (interaction length, and winning process):
235  }
236 
237  // -- [*all processes*] Send operation for processes with well defined XS (see "**!**" above):
238  G4VBiasingOperation* operationToReturn = nullptr;
239  if ( callingProcess->GetWrappedProcess()->GetCurrentInteractionLength() < DBL_MAX/10. ) operationToReturn = fSharedForceInteractionOperation;
240  return operationToReturn;
241 
242 
243  } // -- end of "if ( fCurrentTrackData->fForceCollisionState == ForceCollisionState::toBeForced )"
244 
245 
246  // -- other cases here: particle appearing in the volume by some
247  // -- previous interaction : we decide to not bias these.
248  return 0;
249 
250 }
251 
252 
254  const G4BiasingProcessInterface* /* callingProcess */)
255 {
256  if ( track->GetDefinition() != fParticleToBias ) return nullptr;
257 
258  // -- Apply biasing scheme only to tracks entering the volume.
259  // -- A "cloning" is done:
260  // -- - the primary will be forced flight under a zero weight up to volume exit,
261  // -- where the weight will be restored with proper weight for free flight
262  // -- - the clone will be forced to interact in the volume.
263  if ( track->GetStep()->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ) // -- §§§ extent to case of a track shoot on the boundary ?
264  {
265  // -- check that track is free of undergoing biasing scheme ( no biasing data, or no active active )
266  // -- Get possible track data:
268  if ( fCurrentTrackData != nullptr )
269  {
271  {
272  // -- takes "ownership" (some track data created before, left free, reuse it):
274  }
275  else
276  {
277  // §§§ Would something be really wrong in this case ? Could this be that a process made a zero step ?
278  }
279  }
280  else
281  {
284  }
286  fInitialTrackWeight = track->GetWeight();
288  return fCloningOperation;
289  }
290 
291  // --
292  return nullptr;
293 }
294 
295 
297 {
298  // -- Propose at final state generation the same operation which was proposed at GPIL level,
299  // -- (which is either the force free flight or the force interaction operation).
300  // -- count on the process interface to collect this operation:
301  return callingProcess->GetCurrentOccurenceBiasingOperation();
302 }
303 
304 
306 {
307  fCurrentTrack = track;
308  fCurrentTrackData = nullptr;
309 }
310 
311 
313 {
314  // -- check for consistency, operator should have cleaned the track:
315  if ( fCurrentTrackData != nullptr )
316  {
318  {
320  {
322  ed << "Current track deleted while under biasing by " << GetName() << ". Will result in inconsistencies.";
323  G4Exception(" G4BOptrForceCollision::EndTracking()",
324  "BIAS.GEN.18",
325  JustWarning,
326  ed);
327  }
328  }
329  }
330 }
331 
332 
335  G4VBiasingOperation* operationApplied,
336  const G4VParticleChange* )
337 {
338 
339  if ( fCurrentTrackData == nullptr )
340  {
341  if ( BAC != BAC_None )
342  {
344  ed << " Internal inconsistency : please submit bug report. " << G4endl;
345  G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
346  "BIAS.GEN.20.1",
347  JustWarning,
348  ed);
349  }
350  return;
351  }
352 
354  {
356  auto cloneData = new G4BOptrForceCollisionTrackData( this );
357  cloneData->fForceCollisionState = ForceCollisionState::toBeForced;
359  }
361  {
362  if ( fFreeFlightOperations[callingProcess]->OperationComplete() ) fCurrentTrackData->Reset(); // -- off biasing for this track
363  }
365  {
366  if ( operationApplied != fSharedForceInteractionOperation )
367  {
369  ed << " Internal inconsistency : please submit bug report. " << G4endl;
370  G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
371  "BIAS.GEN.20.2",
372  JustWarning,
373  ed);
374  }
376  {
377  if ( operationApplied != fSharedForceInteractionOperation )
378  {
380  ed << " Internal inconsistency : please submit bug report. " << G4endl;
381  G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
382  "BIAS.GEN.20.3",
383  JustWarning,
384  ed);
385  }
386  }
387  }
388  else
389  {
391  {
393  ed << " Internal inconsistency : please submit bug report. " << G4endl;
394  G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
395  "BIAS.GEN.20.4",
396  JustWarning,
397  ed);
398  }
399  }
400 }
401 
402 
404  G4VBiasingOperation* /*occurenceOperationApplied*/, G4double /*weightForOccurenceInteraction*/,
405  G4VBiasingOperation* finalStateOperationApplied, const G4VParticleChange* /*particleChangeProduced*/ )
406 {
407 
409  {
410  if ( finalStateOperationApplied != fSharedForceInteractionOperation )
411  {
413  ed << " Internal inconsistency : please submit bug report. " << G4endl;
414  G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
415  "BIAS.GEN.20.5",
416  JustWarning,
417  ed);
418  }
419  if ( fSharedForceInteractionOperation->GetInteractionOccured() ) fCurrentTrackData->Reset(); // -- off biasing for this track
420  }
421  else
422  {
424  ed << " Internal inconsistency : please submit bug report. " << G4endl;
425  G4Exception(" G4BOptrForceCollision::OperationApplied(...)",
426  "BIAS.GEN.20.6",
427  JustWarning,
428  ed);
429  }
430 
431 }
432 
G4ParticleDefinition * GetDefinition() const
G4VBiasingOperation * GetCurrentOccurenceBiasingOperation() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4Track * GetCloneTrack() const
void SetAuxiliaryTrackInformation(G4int idx, G4VAuxiliaryTrackInformation *info) const
Definition: G4Track.cc:324
G4String name
Definition: TRTMaterials.hh:40
G4StepStatus GetStepStatus() const
void AddCrossSection(const G4VProcess *, G4double)
const G4BOptrForceCollision * fForceCollisionOperator
G4TrackStatus GetTrackStatus() const
G4BOptrForceCollisionTrackData * fCurrentTrackData
const std::vector< const G4BiasingProcessInterface * > & GetPhysicsBiasingProcessInterfaces() const
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
virtual void EndTracking() final
const G4ThreeVector & GetInitialMomentum() const
const G4Step * GetStep() const
G4ProcessManager * GetProcessManager() const
G4BOptrForceCollision(G4String particleToForce, G4String name="ForceCollision")
G4BiasingAppliedCase
G4VProcess * GetWrappedProcess() const
void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced) final
G4StepPoint * GetPreStepPoint() const
void ResetInitialTrackWeight(G4double w)
G4int GetCurrentStepNumber() const
G4BOptnForceCommonTruncatedExp * fSharedForceInteractionOperation
virtual void StartRun() final
bool G4bool
Definition: G4Types.hh:79
G4double GetCurrentInteractionLength() const
Definition: G4VProcess.hh:462
virtual void StartTracking(const G4Track *track) final
G4BOptnCloning * fCloningOperation
virtual void Configure() final
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
G4VAuxiliaryTrackInformation * GetAuxiliaryTrackInformation(G4int idx) const
Definition: G4Track.cc:340
std::map< const G4BiasingProcessInterface *, G4BOptnForceFreeFlight * > fFreeFlightOperations
void SetCloneWeights(G4double clone1Weight, G4double clone2Weight)
G4ThreeVector GetMomentum() const
static G4ParticleTable * GetParticleTable()
static G4int Register(const G4String &)
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
#define DBL_MIN
Definition: templates.hh:75
G4double GetWeight() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
#define DBL_MAX
Definition: templates.hh:83
virtual void ConfigureForWorker() final
const G4String GetName() const
const G4ParticleDefinition * fParticleToBias
const G4BiasingProcessSharedData * GetSharedData() const