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