Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BiasingProcessInterface.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// G4BiasingProcessInterface
27// --------------------------------------------------------------------
28
30#include "G4VBiasingOperator.hh"
36#include "G4ProcessManager.hh"
39
40G4Cache<G4bool> G4BiasingProcessInterface::fResetInteractionLaws;// = true;
41G4Cache<G4bool> G4BiasingProcessInterface::fCommonStart;// = true;
42G4Cache<G4bool> G4BiasingProcessInterface::fCommonEnd;// = true;
43G4Cache<G4bool> G4BiasingProcessInterface::fDoCommonConfigure;
44
46 : G4VProcess( name ),
47 fResetWrappedProcessInteractionLength( true )
48{
49 for (G4int i = 0 ; i < 8 ; i++) fFirstLastFlags[i] = false;
50 fResetInteractionLaws.Put( true );
51 fCommonStart .Put( true );
52 fCommonEnd .Put( true );
53 fDoCommonConfigure .Put( true );
54}
55
56
59 G4bool wrappedIsAtRest,
60 G4bool wrappedIsAlongStep,
61 G4bool wrappedIsPostStep,
62 G4String useThisName)
63 : G4VProcess(useThisName != ""
64 ? useThisName
65 : "biasWrapper("+wrappedProcess->GetProcessName()+")",
66 wrappedProcess->GetProcessType()),
67 fWrappedProcess( wrappedProcess ),
68 fIsPhysicsBasedBiasing( true ),
69 fWrappedProcessIsAtRest( wrappedIsAtRest ),
70 fWrappedProcessIsAlong( wrappedIsAlongStep ),
71 fWrappedProcessIsPost( wrappedIsPostStep )
72{
73 for (G4int i = 0 ; i < 8 ; ++i)
74 fFirstLastFlags[i] = false;
75 fResetInteractionLaws.Put( true );
76 fCommonStart.Put(true);
77 fCommonEnd.Put(true);
78 fDoCommonConfigure.Put(true);
79
80 SetProcessSubType(fWrappedProcess->GetProcessSubType());
81
82 // -- create physical interaction law:
83 fPhysicalInteractionLaw = new G4InteractionLawPhysical("PhysicalInteractionLawFor("+GetProcessName()+")");
84 // -- instantiate particle change wrapper for occurrence biaising:
85 fOccurenceBiasingParticleChange = new G4ParticleChangeForOccurenceBiasing("biasingPCfor"+GetProcessName());
86 // -- instantiate a "do nothing" particle change:
87 fDummyParticleChange = new G4ParticleChangeForNothing();
88}
89
91{
92 delete fPhysicalInteractionLaw;
93 delete fOccurenceBiasingParticleChange;
94 delete fDummyParticleChange;
95}
96
99{
100 const auto & itr = G4BiasingProcessSharedData::fSharedDataMap.Find( mgr );
101 if ( itr != G4BiasingProcessSharedData::fSharedDataMap.End( ) )
102 {
103 return (*itr).second;
104 }
105 else return nullptr;
106}
107
108
110{
111 fCurrentTrack = track;
112 if ( fIsPhysicsBasedBiasing ) fWrappedProcess->StartTracking(fCurrentTrack);
113 fOccurenceBiasingOperation = nullptr;
114 fPreviousOccurenceBiasingOperation = nullptr;
115 fFinalStateBiasingOperation = nullptr;
116 fPreviousFinalStateBiasingOperation = nullptr;
117 fNonPhysicsBiasingOperation = nullptr;
118 fPreviousNonPhysicsBiasingOperation = nullptr;
119 fBiasingInteractionLaw = nullptr;
120 fPreviousBiasingInteractionLaw = nullptr;
121
122 fPreviousStepSize = -1.0;
123
124 fResetWrappedProcessInteractionLength = false;
125
126 if ( fCommonStart.Get() )
127 {
128 fCommonStart.Put( false );// = false;
129 fCommonEnd.Put( true );// = true;
130
131 fSharedData->fCurrentBiasingOperator = nullptr;
132 fSharedData->fPreviousBiasingOperator = nullptr;
133
134 // -- Add a "fSharedData->nStarting" here and outside bracket "fSharedData->nStarting++" and " if (fSharedData->nStarting) == fSharedData->(vector interface length)"
135 // -- call to the loop "StartTracking" of operators"
136
137 for (std::size_t optr=0 ; optr<(G4VBiasingOperator::GetBiasingOperators()).size(); ++optr)
138 {
139 (G4VBiasingOperator::GetBiasingOperators())[optr]->StartTracking( fCurrentTrack );
140 }
141 }
142}
143
145{
146 if ( fIsPhysicsBasedBiasing )
147 fWrappedProcess->EndTracking();
148 if ( fSharedData->fCurrentBiasingOperator)
149 (fSharedData->fCurrentBiasingOperator)->ExitingBiasing(fCurrentTrack, this);
150 fBiasingInteractionLaw = nullptr;
151
152 // -- Inform operators of end of tracking:
153 if ( fCommonEnd.Get() )
154 {
155 fCommonEnd .Put( false );// = false;
156 fCommonStart.Put( true );// = true;
157
158 for ( std::size_t optr=0; optr<(G4VBiasingOperator::GetBiasingOperators()).size(); ++optr)
159 {
161 }
162 // -- for above loop, do as in StartTracking.
163 }
164}
165
168 G4double previousStepSize,
170{
171 // ---------------------------------------------------------------------------------------------------
172 // -- The "biasing process master" takes care of updating the biasing operator, and for all biasing
173 // -- processes it invokes the PostStepGPIL of physical wrapped processes (anticipate stepping manager
174 // -- call ! ) to make all cross-sections updated with current step, and hence available before the
175 // -- first call to the biasing operator.
176 // ---------------------------------------------------------------------------------------------------
177 if ( fIamFirstGPIL )
178 {
179 // -- Update previous biasing operator, and assume the operator stays the same by
180 // -- default and that it is not left at the beginning of this step. These
181 // -- assumptions might be wrong if there is a volume change (in paralllel or
182 // -- mass geometries) in what case the flags will be updated.
183 fSharedData->fPreviousBiasingOperator = fSharedData->fCurrentBiasingOperator;
184 fSharedData->fIsNewOperator = false;
185 fSharedData->fLeavingPreviousOperator = false;
186 // -- If new volume, either in mass or parallel geometries, get possible new biasing operator:
187 // -------------------------------------------------------------------------------------------
188 // -- Get biasing operator in parallel geometries:
189 G4bool firstStepInParallelVolume = false;
190 if ( fSharedData->fParallelGeometriesLimiterProcess )
191 {
192 G4VBiasingOperator* newParallelOperator( nullptr );
193 G4bool firstStep = ( track.GetCurrentStepNumber() == 1 );
194 std::size_t iParallel = 0;
195 for ( auto wasLimiting : fSharedData->fParallelGeometriesLimiterProcess->GetWasLimiting() )
196 {
197 if ( firstStep || wasLimiting )
198 {
199 firstStepInParallelVolume = true;
200 auto tmpParallelOperator = G4VBiasingOperator::
201 GetBiasingOperator((fSharedData->fParallelGeometriesLimiterProcess
202 ->GetCurrentVolumes()[iParallel])
203 ->GetLogicalVolume());
204 if ( newParallelOperator )
205 {
206 if ( tmpParallelOperator )
207 {
209 ed << " Several biasing operators are defined at the same place\n"
210 << " in parallel geometries ! Found:\n";
211 ed << " - `" << newParallelOperator->GetName() << "' and \n";
212 ed << " - `" << tmpParallelOperator->GetName() << "'.\n";
213 ed << " Keeping `" << newParallelOperator->GetName()
214 << "'. Behavior not guaranteed ! Please consider having only one operator at a place."
215 << G4endl;
216 G4Exception(" G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(...)",
217 "BIAS.GEN.30", JustWarning, ed);
218 }
219 }
220 else newParallelOperator = tmpParallelOperator;
221 }
222 ++iParallel;
223 }
224 fSharedData->fParallelGeometryOperator = newParallelOperator;
225 } // -- end of " if ( fSharedData->fParallelGeometriesLimiterProcess )"
226
227 // -- Get biasing operator in mass geometry:
228 // -- [ยงยง Note : bug with this first step ? Does not work if previous step was concurrently limited with geometry. Might make use of safety at last point ?]
229 G4bool firstStepInVolume = ( (track.GetStep()->GetPreStepPoint()->GetStepStatus() == fGeomBoundary)
230 || (track.GetCurrentStepNumber() == 1) );
231 // fSharedData->fIsNewOperator = false;
232 // fSharedData->fLeavingPreviousOperator = false;
233 if ( firstStepInVolume )
234 {
237 fSharedData->fMassGeometryOperator = newOperator;
238 if ( ( newOperator != nullptr ) && ( fSharedData->fParallelGeometryOperator != nullptr ) )
239 {
241 ed << " Biasing operators are defined at the same place in mass and parallel geometries ! Found:\n";
242 ed << " - `" << fSharedData->fParallelGeometryOperator->GetName() << "' in parallel geometry and \n";
243 ed << " - `" << newOperator->GetName() << "' in mass geometry.\n";
244 ed << " Keeping `" << fSharedData->fParallelGeometryOperator->GetName() << "'. Behavior not guaranteed ! Please consider having only one operator at a place. " << G4endl;
245 G4Exception(" G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(...)",
246 "BIAS.GEN.31", JustWarning, ed);
247 }
248 }
249
250 // -- conclude the operator selection, giving priority to parallel geometry (as told in exception message BIAS.GEN.30):
251 if ( firstStepInVolume || firstStepInParallelVolume )
252 {
253 G4VBiasingOperator* newOperator = fSharedData->fParallelGeometryOperator;
254 if ( newOperator == nullptr )
255 newOperator = fSharedData->fMassGeometryOperator;
256
257 fSharedData->fCurrentBiasingOperator = newOperator ;
258
259 if ( newOperator != fSharedData->fPreviousBiasingOperator )
260 {
261 fSharedData->fLeavingPreviousOperator = ( fSharedData->fPreviousBiasingOperator != nullptr ) ;
262 fSharedData->fIsNewOperator = ( newOperator != nullptr );
263 }
264 }
265
266 // -- calls to wrapped process PostStepGPIL's:
267 // -------------------------------------------
268 // -- Each physics wrapper process has its
269 // -- fWrappedProcessPostStepGPIL ,
270 // -- fWrappedProcessForceCondition ,
271 // -- fWrappedProcessInteractionLength
272 // -- updated.
273 if ( fSharedData->fCurrentBiasingOperator != nullptr )
274 {
275 for (std::size_t i=0; i<(fSharedData->fPhysicsBiasingProcessInterfaces).size(); ++i)
276 {
277 (fSharedData->fPhysicsBiasingProcessInterfaces)[i]->InvokeWrappedProcessPostStepGPIL( track, previousStepSize, condition );
278 }
279 }
280 } // -- end of "if ( fIamFirstGPIL )"
281
282 // -- Remember previous operator and proposed operations, if any, and reset:
283 // -------------------------------------------------------------------------
284 // -- remember only in case some biasing might be called
285 if ( ( fSharedData->fPreviousBiasingOperator != nullptr ) ||
286 ( fSharedData->fCurrentBiasingOperator != nullptr ) )
287 {
288 fPreviousOccurenceBiasingOperation = fOccurenceBiasingOperation;
289 fPreviousFinalStateBiasingOperation = fFinalStateBiasingOperation;
290 fPreviousNonPhysicsBiasingOperation = fNonPhysicsBiasingOperation;
291 fPreviousBiasingInteractionLaw = fBiasingInteractionLaw;
292 // -- reset:
293 fOccurenceBiasingOperation = nullptr;
294 fFinalStateBiasingOperation = nullptr;
295 fNonPhysicsBiasingOperation = nullptr;
296 fBiasingInteractionLaw = nullptr;
297 // -- Physics PostStep and AlongStep GPIL
298 // fWrappedProcessPostStepGPIL : updated by InvokeWrappedProcessPostStepGPIL(...) above
299 fBiasingPostStepGPIL = DBL_MAX;
300 // fWrappedProcessInteractionLength : updated by InvokeWrappedProcessPostStepGPIL(...) above; inverse of analog cross-section.
301 // fWrappedProcessForceCondition : updated by InvokeWrappedProcessPostStepGPIL(...) above
302 fBiasingForceCondition = NotForced;
303 fWrappedProcessAlongStepGPIL = DBL_MAX;
304 fBiasingAlongStepGPIL = DBL_MAX;
305 fWrappedProcessGPILSelection = NotCandidateForSelection;
306 fBiasingGPILSelection = NotCandidateForSelection;
307 // -- for helper:
308 fPreviousStepSize = previousStepSize;
309 }
310
311 // -- previous step size value; it is switched to zero if resetting a wrapped process:
312 // -- (same trick used than in InvokedWrappedProcessPostStepGPIL )
313 G4double usedPreviousStepSize = previousStepSize;
314
315 // ----------------------------------------------
316 // -- If leaving a biasing operator, let it know:
317 // ----------------------------------------------
318 if ( fSharedData->fLeavingPreviousOperator )
319 {
320 (fSharedData->fPreviousBiasingOperator)->ExitingBiasing( &track, this );
321 // -- if no further biasing operator, reset process behavior to standard tracking:
322 if ( fSharedData->fCurrentBiasingOperator == nullptr )
323 {
324 ResetForUnbiasedTracking();
325 if ( fIsPhysicsBasedBiasing )
326 {
327 // -- if the physics process has been under occurrence biasing, reset it:
328 if ( fResetWrappedProcessInteractionLength )
329 {
330 fResetWrappedProcessInteractionLength = false;
331 fWrappedProcess->ResetNumberOfInteractionLengthLeft();
332 // -- We set "previous step size" as 0.0, to let the process believe this is first step:
333 usedPreviousStepSize = 0.0;
334 }
335 }
336 }
337 }
338
339 // --------------------------------------------------------------
340 // -- no operator : analog tracking if physics-based, or nothing:
341 // --------------------------------------------------------------
342 if ( fSharedData->fCurrentBiasingOperator == nullptr )
343 {
344 // -- take note of the "usedPreviousStepSize" value:
345 if ( fIsPhysicsBasedBiasing )
346 {
347 return fWrappedProcess->PostStepGetPhysicalInteractionLength(track, usedPreviousStepSize, condition);
348 }
349 else
350 {
352 return DBL_MAX;
353 }
354 }
355
356 // --------------------------------------------------
357 // -- A biasing operator exists. Proceed with
358 // -- treating non-physics and physics biasing cases:
359 //---------------------------------------------------
360
361 // -- non-physics-based biasing case:
362 // ----------------------------------
363 if ( !fIsPhysicsBasedBiasing )
364 {
365 fNonPhysicsBiasingOperation = (fSharedData->fCurrentBiasingOperator)->GetProposedNonPhysicsBiasingOperation( &track, this );
366 if ( fNonPhysicsBiasingOperation == nullptr )
367 {
369 return DBL_MAX;
370 }
371 return fNonPhysicsBiasingOperation->DistanceToApplyOperation(&track, previousStepSize, condition);
372 }
373
374 // -- Physics based biasing case:
375 // ------------------------------
376 // -- Ask for possible GPIL biasing operation:
377 fOccurenceBiasingOperation = (fSharedData->fCurrentBiasingOperator)->GetProposedOccurenceBiasingOperation( &track, this );
378
379 // -- no operation for occurrence biasing, analog GPIL returns the wrapped process GPIL and condition values
380 if ( fOccurenceBiasingOperation == nullptr )
381 {
382 *condition = fWrappedProcessForceCondition;
383 return fWrappedProcessPostStepGPIL;
384 }
385
386 // -- A valid GPIL biasing operation has been proposed:
387 // -- 0) remember wrapped process will need to be reset on biasing exit, if particle survives:
388 fResetWrappedProcessInteractionLength = true;
389 // -- 1) update process interaction length for reference analog interaction law ( fWrappedProcessInteractionLength updated/collected above):
390 fPhysicalInteractionLaw->SetPhysicalCrossSection( 1.0 / fWrappedProcessInteractionLength );
391 // -- 2) Collect biasing interaction law:
392 // -- The interaction law pointer is collected as a const pointer to the interaction law object.
393 // -- This interaction law will be kept under control of the biasing operation, which is the only
394 // -- entity that will change the state of the biasing interaction law.
395 // -- The force condition for biasing is asked at the same time, passing the analog one as default:
396 fBiasingForceCondition = fWrappedProcessForceCondition;
397 fBiasingInteractionLaw = fOccurenceBiasingOperation->ProvideOccurenceBiasingInteractionLaw( this, fBiasingForceCondition );
398 // -- 3) Ask operation to sample the biasing interaction law:
399 fBiasingPostStepGPIL = fBiasingInteractionLaw->GetSampledInteractionLength();
400
401 // -- finish
402 *condition = fBiasingForceCondition;
403 return fBiasingPostStepGPIL;
404}
405
407 const G4Step& step)
408{
409 // ---------------------------------------
410 // -- case outside of volume with biasing:
411 // ---------------------------------------
412 if ( fSharedData->fCurrentBiasingOperator == nullptr )
413 return fWrappedProcess->PostStepDoIt(track, step);
414
415 // ----------------------------
416 // -- non-physics biasing case:
417 // ----------------------------
418 if ( !fIsPhysicsBasedBiasing )
419 {
420 G4VParticleChange* particleChange = fNonPhysicsBiasingOperation->GenerateBiasingFinalState( &track, &step );
421 (fSharedData->fCurrentBiasingOperator)->ReportOperationApplied( this, BAC_NonPhysics, fNonPhysicsBiasingOperation, particleChange );
422 return particleChange;
423 }
424
425 // -- physics biasing case:
426 // ------------------------
427 // -- It proceeds with the following logic:
428 // -- 1) Obtain the final state
429 // -- This final state may be analog or biased.
430 // -- The biased final state is obtained through a biasing operator
431 // -- returned by the operator.
432 // -- 2) The biased final state may be asked to be "force as it is"
433 // -- in what case the particle change is returned as is to the
434 // -- stepping.
435 // -- In all other cases (analog final state or biased final but
436 // -- not forced) the final state weight may be modified by the
437 // -- occurrence biasing, if such an occurrence biasing is at play.
438
439 // -- Get final state, biased or analog:
440 G4VParticleChange* finalStateParticleChange;
442 fFinalStateBiasingOperation = (fSharedData->fCurrentBiasingOperator)->GetProposedFinalStateBiasingOperation( &track, this );
443 // -- Flag below is to force the biased generated particle change to be returned "as is" to the stepping, disregarding there
444 // -- was or not a occurrence biasing that would apply. Weight relevance under full responsibility of the biasing operation.
445 G4bool forceBiasedFinalState = false;
446 if ( fFinalStateBiasingOperation != nullptr )
447 {
448 finalStateParticleChange = fFinalStateBiasingOperation->ApplyFinalStateBiasing( this, &track, &step, forceBiasedFinalState );
449 BAC = BAC_FinalState;
450 }
451 else
452 {
453 finalStateParticleChange = fWrappedProcess->PostStepDoIt(track, step);
454 BAC = BAC_None ;
455 }
456
457 // -- if no occurrence biasing operation, we're done:
458 if ( fOccurenceBiasingOperation == nullptr )
459 {
460 (fSharedData->fCurrentBiasingOperator)->ReportOperationApplied( this, BAC, fFinalStateBiasingOperation, finalStateParticleChange );
461 return finalStateParticleChange;
462 }
463
464 // -- if biased final state has been asked to be forced, we're done:
465 if ( forceBiasedFinalState )
466 {
467 (fSharedData->fCurrentBiasingOperator)->ReportOperationApplied( this, BAC, fFinalStateBiasingOperation, finalStateParticleChange );
468 return finalStateParticleChange;
469 }
470
471 // -- If occurrence biasing, applies the occurrence biasing weight correction on top of final state (biased or not):
472 G4double weightForInteraction = 1.0;
473 if ( !fBiasingInteractionLaw->IsSingular() )
474 {
475 weightForInteraction = fPhysicalInteractionLaw->ComputeEffectiveCrossSectionAt(step.GetStepLength())
476 / fBiasingInteractionLaw ->ComputeEffectiveCrossSectionAt(step.GetStepLength());
477 }
478 else
479 {
480 // -- at this point effective XS can only be infinite, if not, there is a logic problem
481 if ( !fBiasingInteractionLaw->IsEffectiveCrossSectionInfinite() )
482 {
484 ed << "Internal inconsistency in cross-section handling. Please report !" << G4endl;
485 G4Exception(" G4BiasingProcessInterface::PostStepDoIt(...)",
486 "BIAS.GEN.02", JustWarning, ed);
487 // -- if XS is infinite, weight is zero (and will stay zero), but we'll do differently.
488 // -- Should foresee in addition something to remember that in case of singular
489 // -- distribution, weight can only be partly calculated
490 }
491 }
492
493 if ( weightForInteraction <= 0. )
494 {
496 ed << " Negative interaction weight : w_I = "
497 << weightForInteraction << " XS_I(phys) = "
498 << fBiasingInteractionLaw ->ComputeEffectiveCrossSectionAt(step.GetStepLength())
499 <<" XS_I(bias) = "
500 << fPhysicalInteractionLaw->ComputeEffectiveCrossSectionAt(step.GetStepLength())
501 << " step length = " << step.GetStepLength()
502 << " Interaction law = `" << fBiasingInteractionLaw << "'"
503 << G4endl;
504 G4Exception(" G4BiasingProcessInterface::PostStepDoIt(...)",
505 "BIAS.GEN.03", JustWarning, ed);
506 }
507
508 (fSharedData->fCurrentBiasingOperator)
509 ->ReportOperationApplied( this, BAC, fOccurenceBiasingOperation,
510 weightForInteraction,
511 fFinalStateBiasingOperation,
512 finalStateParticleChange );
513
514 fOccurenceBiasingParticleChange->SetOccurenceWeightForInteraction( weightForInteraction );
515 fOccurenceBiasingParticleChange->SetSecondaryWeightByProcess( true );
516 fOccurenceBiasingParticleChange->SetWrappedParticleChange( finalStateParticleChange );
517 fOccurenceBiasingParticleChange->ProposeTrackStatus( finalStateParticleChange->GetTrackStatus() );
518 fOccurenceBiasingParticleChange->StealSecondaries(); // -- this also makes weightForInteraction applied to secondaries stolen
519
520 // -- finish:
521 return fOccurenceBiasingParticleChange;
522}
523
524// -- AlongStep methods:
527 G4double previousStepSize,
528 G4double currentMinimumStep,
529 G4double& proposedSafety,
530 G4GPILSelection* selection)
531{
532 // -- for helper methods:
533 fCurrentMinimumStep = currentMinimumStep;
534 fProposedSafety = proposedSafety;
535
536 // -- initialization default case:
537 fWrappedProcessAlongStepGPIL = DBL_MAX;
538 *selection = NotCandidateForSelection;
539 // ---------------------------------------
540 // -- case outside of volume with biasing:
541 // ---------------------------------------
542 if ( fSharedData->fCurrentBiasingOperator == nullptr )
543 {
544 if ( fWrappedProcessIsAlong )
545 fWrappedProcessAlongStepGPIL = fWrappedProcess
546 ->AlongStepGetPhysicalInteractionLength(track, previousStepSize,
547 currentMinimumStep,
548 proposedSafety, selection);
549 return fWrappedProcessAlongStepGPIL;
550 }
551
552 // --------------------------------------------------------------------
553 // -- non-physics based biasing: no along operation expected (for now):
554 // --------------------------------------------------------------------
555 if ( !fIsPhysicsBasedBiasing ) return fWrappedProcessAlongStepGPIL;
556
557 // ----------------------
558 // -- physics-based case:
559 // ----------------------
560 if ( fOccurenceBiasingOperation == nullptr )
561 {
562 if ( fWrappedProcessIsAlong )
563 fWrappedProcessAlongStepGPIL = fWrappedProcess
564 ->AlongStepGetPhysicalInteractionLength(track, previousStepSize,
565 currentMinimumStep,
566 proposedSafety, selection);
567 return fWrappedProcessAlongStepGPIL;
568 }
569
570 // -----------------------------------------------------------
571 // -- From here we have an valid occurrence biasing operation:
572 // -----------------------------------------------------------
573 // -- Give operation opportunity to shorten step proposed by physics process:
574 fBiasingAlongStepGPIL = fOccurenceBiasingOperation->ProposeAlongStepLimit( this );
575 G4double minimumStep = fBiasingAlongStepGPIL < currentMinimumStep
576 ? fBiasingAlongStepGPIL : currentMinimumStep;
577 // -- wrapped process is called with minimum step ( <= currentMinimumStep passed ) : an along process can not
578 // -- have its operation stretched over what it expects:
579 if ( fWrappedProcessIsAlong )
580 {
581 fWrappedProcessAlongStepGPIL = fWrappedProcess
582 ->AlongStepGetPhysicalInteractionLength(track, previousStepSize,
583 minimumStep,
584 proposedSafety, selection);
585 fWrappedProcessGPILSelection = *selection;
586 fBiasingGPILSelection = fOccurenceBiasingOperation
587 ->ProposeGPILSelection( fWrappedProcessGPILSelection );
588 }
589 else
590 {
591 fBiasingGPILSelection = fOccurenceBiasingOperation
592 ->ProposeGPILSelection( NotCandidateForSelection );
593 fWrappedProcessAlongStepGPIL = fBiasingAlongStepGPIL;
594 }
595
596 *selection = fBiasingGPILSelection;
597
598 return fWrappedProcessAlongStepGPIL;
599}
600
603 const G4Step& step)
604{
605 // ---------------------------------------
606 // -- case outside of volume with biasing:
607 // ---------------------------------------
608 if ( fSharedData->fCurrentBiasingOperator == nullptr )
609 {
610 if ( fWrappedProcessIsAlong )
611 {
612 return fWrappedProcess->AlongStepDoIt(track, step);
613 }
614 else
615 {
616 fDummyParticleChange->Initialize( track );
617 return fDummyParticleChange;
618 }
619 }
620
621 // -----------------------------------
622 // -- case inside volume with biasing:
623 // -----------------------------------
624 if ( fWrappedProcessIsAlong )
625 {
626 fOccurenceBiasingParticleChange
627 ->SetWrappedParticleChange(fWrappedProcess->AlongStepDoIt(track, step));
628 }
629 else
630 {
631 fOccurenceBiasingParticleChange->SetWrappedParticleChange ( nullptr );
632 fOccurenceBiasingParticleChange->ProposeTrackStatus( track.GetTrackStatus() );
633 }
634 G4double weightForNonInteraction (1.0);
635 if ( fBiasingInteractionLaw != nullptr )
636 {
637 weightForNonInteraction =
638 fPhysicalInteractionLaw->ComputeNonInteractionProbabilityAt(step.GetStepLength()) /
639 fBiasingInteractionLaw ->ComputeNonInteractionProbabilityAt(step.GetStepLength());
640
641 fOccurenceBiasingOperation->AlongMoveBy( this, &step, weightForNonInteraction );
642
643 if ( weightForNonInteraction <= 0. )
644 {
646 ed << " Negative non interaction weight : w_NI = " << weightForNonInteraction <<
647 " p_NI(phys) = " << fPhysicalInteractionLaw->ComputeNonInteractionProbabilityAt(step.GetStepLength()) <<
648 " p_NI(bias) = " << fBiasingInteractionLaw ->ComputeNonInteractionProbabilityAt(step.GetStepLength()) <<
649 " step length = " << step.GetStepLength() <<
650 " biasing interaction law = `" << fBiasingInteractionLaw->GetName() << "'" << G4endl;
651 G4Exception(" G4BiasingProcessInterface::AlongStepDoIt(...)",
652 "BIAS.GEN.04", JustWarning, ed);
653 }
654 }
655
656 fOccurenceBiasingParticleChange
657 ->SetOccurenceWeightForNonInteraction( weightForNonInteraction );
658
659 return fOccurenceBiasingParticleChange;
660}
661
662// -- AtRest methods
666{
667 return fWrappedProcess->AtRestGetPhysicalInteractionLength(track, condition);
668}
669
671 const G4Step& step)
672{
673 return fWrappedProcess->AtRestDoIt(track, step);
674}
675
677{
678 if ( fWrappedProcess != nullptr )
679 return fWrappedProcess->IsApplicable(pd);
680 else
681 return true;
682}
683
685{
686 // -- Master for this process:
688 // -- Master for wrapped process:
689 if ( fWrappedProcess != nullptr )
690 {
691 const G4BiasingProcessInterface* thisWrapperMaster
693 // -- paranoia check: (?)
694 G4VProcess* wrappedMaster = nullptr;
695 wrappedMaster = thisWrapperMaster->GetWrappedProcess();
696 fWrappedProcess->SetMasterProcess( wrappedMaster );
697 }
698}
699
702{
703 // -- Sequential mode : called second (after PreparePhysicsTable(..))
704 // -- MT mode : called second (after PreparePhysicsTable(..)) by master thread.
705 // -- Corresponding process instance not used then by tracking.
706 // -- PreparePhysicsTable(...) has been called first for all processes,
707 // -- so the first/last flags and G4BiasingProcessInterface vector of processes have
708 // -- been properly setup, fIamFirstGPIL is valid.
709 if ( fWrappedProcess != nullptr )
710 {
711 fWrappedProcess->BuildPhysicsTable(pd);
712 }
713
714 if ( fIamFirstGPIL )
715 {
716 // -- Re-order vector of processes to match that of the GPIL
717 // -- (made for fIamFirstGPIL, but important is to have it made once):
718 ReorderBiasingVectorAsGPIL();
719 // -- Let operators to configure themselves for the master thread or for sequential mode.
720 // -- Intended here is in particular the registration to physics model catalog.
721 // -- The fDoCommonConfigure is to ensure that this Configure is made by only one process (othewise each first process makes the call):
722 if ( fDoCommonConfigure.Get() )
723 {
724 for ( std::size_t optr=0; optr<(G4VBiasingOperator::GetBiasingOperators()).size(); ++optr)
725 {
726 (G4VBiasingOperator::GetBiasingOperators())[optr]->Configure( );
727 }
728 fDoCommonConfigure.Put(false);
729 }
730 }
731}
732
735{
736 // -- Sequential mode : called first (before BuildPhysicsTable(..))
737 // -- MT mode : called first (before BuildPhysicsTable(..)) by master thread.
738 // -- Corresponding process instance not used then by tracking.
739 // -- Let process finding its first/last position in the process manager:
740 SetUpFirstLastFlags();
741 if ( fWrappedProcess != nullptr )
742 {
743 fWrappedProcess->PreparePhysicsTable(pd);
744 }
745}
746
749{
750 if ( fWrappedProcess != nullptr )
751 return fWrappedProcess->StorePhysicsTable(pd, s, f);
752 else
753 return false;
754}
755
758{
759 if ( fWrappedProcess != nullptr )
760 return fWrappedProcess->RetrievePhysicsTable(pd, s, f);
761 else
762 return false;
763}
764
766{
767 if ( fWrappedProcess != nullptr )
768 fWrappedProcess->SetProcessManager(mgr);
769 else
771
772 // -- initialize fSharedData pointer:
773 if (G4BiasingProcessSharedData::fSharedDataMap.Find(mgr)
774 == G4BiasingProcessSharedData::fSharedDataMap.End() )
775 {
776 fSharedData = new G4BiasingProcessSharedData( mgr );
777 G4BiasingProcessSharedData::fSharedDataMap[mgr] = fSharedData;
778 }
779 else
780 {
781 fSharedData = G4BiasingProcessSharedData::fSharedDataMap[mgr] ;
782 }
783 // -- augment list of co-operating processes:
784 fSharedData->fBiasingProcessInterfaces.push_back( this );
785 fSharedData->fPublicBiasingProcessInterfaces.push_back( this );
786 if ( fIsPhysicsBasedBiasing )
787 {
788 fSharedData->fPhysicsBiasingProcessInterfaces.push_back( this );
789 fSharedData-> fPublicPhysicsBiasingProcessInterfaces.push_back( this );
790 }
791 else
792 {
793 fSharedData->fNonPhysicsBiasingProcessInterfaces.push_back( this );
794 fSharedData->fPublicNonPhysicsBiasingProcessInterfaces.push_back( this );
795 }
796 // -- remember process manager:
797 fProcessManager = mgr;
798}
799
801{
802 if ( fWrappedProcess != nullptr )
803 return fWrappedProcess->GetProcessManager();
804 else
806}
807
810{
811 // -- Sequential mode : not called
812 // -- MT mode : called after PrepareWorkerPhysicsTable(..)
813 // -- PrepareWorkerPhysicsTable(...) has been called first for all processes,
814 // -- so the first/last flags and G4BiasingProcessInterface vector of processes have
815 // -- been properly setup, fIamFirstGPIL is valid.
816 if ( fWrappedProcess != nullptr )
817 {
818 fWrappedProcess->BuildWorkerPhysicsTable(pd);
819 }
820
821 if ( fIamFirstGPIL )
822 {
823 // -- Re-order vector of processes to match that of the GPIL
824 // -- (made for fIamFirstGPIL, but important is to have it made once):
825 ReorderBiasingVectorAsGPIL();
826 // -- Let operators to configure themselves for the worker thread, if needed.
827 // -- Registration to physics model catalog **IS NOT** to be made here, but in Configure().
828 // -- The fDoCommonConfigure is to ensure that this Configure is made by only one process (othewise each first process makes the call):
829 if ( fDoCommonConfigure.Get() )
830 {
831 for ( std::size_t optr=0 ; optr<(G4VBiasingOperator::GetBiasingOperators()).size(); ++optr)
832 {
833 (G4VBiasingOperator::GetBiasingOperators())[optr]->ConfigureForWorker( );
834 }
835 fDoCommonConfigure.Put(false);
836 }
837 }
838}
839
842{
843 // -- Sequential mode : not called
844 // -- MT mode : called first, before BuildWorkerPhysicsTable(..)
845 // -- Let process finding its first/last position in the process manager:
846 SetUpFirstLastFlags();
847
848 if ( fWrappedProcess != nullptr )
849 {
850 fWrappedProcess->PrepareWorkerPhysicsTable(pd);
851 }
852}
853
855{
856 if ( fWrappedProcess != nullptr )
857 fWrappedProcess->ResetNumberOfInteractionLengthLeft();
858}
859
862{
863 G4int iPhys = ( physOnly ) ? 1 : 0;
864 return fFirstLastFlags[IdxFirstLast( 1, 1, iPhys)];
865}
866
869{
870 G4int iPhys = ( physOnly ) ? 1 : 0;
871 return fFirstLastFlags[IdxFirstLast( 0, 1, iPhys)];
872}
873
876{
877 G4int iPhys = ( physOnly ) ? 1 : 0;
878 return fFirstLastFlags[IdxFirstLast( 1, 0, iPhys)];
879}
880
883{
884 G4int iPhys = ( physOnly ) ? 1 : 0;
885 return fFirstLastFlags[IdxFirstLast( 0, 0, iPhys)];
886}
887
890{
891 G4bool isFirst = true;
892 const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeGPIL);
893 G4int thisIdx(-1);
894 for ( auto i = 0; i < (G4int)pv->size(); ++i )
895 {
896 if ( (*pv)(i) == this ) { thisIdx = i; break; }
897 }
898 if ( thisIdx < 0 ) return false; // -- to ignore pure along processes
899 for ( std::size_t i=0; i<(fSharedData->fBiasingProcessInterfaces).size(); ++i )
900 {
901 if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
902 {
903 G4int thatIdx(-1);
904 for (auto j = 0; j < (G4int)pv->size(); ++j )
905 {
906 if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] )
907 {
908 thatIdx = j; break;
909 }
910 }
911 if ( thatIdx >= 0 ) // -- to ignore pure along processes
912 {
913 if ( thisIdx > thatIdx )
914 {
915 isFirst = false;
916 break;
917 }
918 }
919 }
920 }
921 return isFirst;
922}
923
926{
927 G4bool isLast = true;
928 const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeGPIL);
929 G4int thisIdx(-1);
930 for (auto i = 0; i < (G4int)pv->size(); ++i )
931 {
932 if ( (*pv)(i) == this ) { thisIdx = i; break; }
933 }
934 if ( thisIdx < 0 ) return false; // -- to ignore pure along processes
935 for (std::size_t i=0; i<(fSharedData->fBiasingProcessInterfaces).size(); ++i)
936 {
937 if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
938 {
939 G4int thatIdx(-1);
940 for (auto j = 0; j < (G4int)pv->size(); ++j )
941 {
942 if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] )
943 {
944 thatIdx = j; break;
945 }
946 }
947 if ( thatIdx >= 0 ) // -- to ignore pure along processes
948 {
949 if ( thisIdx < thatIdx )
950 {
951 isLast = false;
952 break;
953 }
954 }
955 }
956 }
957 return isLast;
958}
959
962{
963 G4bool isFirst = true;
964 const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeDoIt);
965 G4int thisIdx(-1);
966 for (auto i = 0; i < (G4int)pv->size(); ++i )
967 {
968 if ( (*pv)(i) == this ) { thisIdx = i; break; }
969 }
970 if ( thisIdx < 0 ) return false; // -- to ignore pure along processes
971 for (std::size_t i=0; i<(fSharedData->fBiasingProcessInterfaces).size(); ++i)
972 {
973 if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
974 {
975 G4int thatIdx(-1);
976 for (auto j = 0; j < (G4int)pv->size(); ++j )
977 {
978 if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] )
979 {
980 thatIdx = j; break;
981 }
982 }
983 if ( thatIdx >= 0 ) // -- to ignore pure along processes
984 {
985 if ( thisIdx > thatIdx )
986 {
987 isFirst = false;
988 break;
989 }
990 }
991 }
992 }
993 return isFirst;
994}
995
998{
999 G4bool isLast = true;
1000 const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeDoIt);
1001 G4int thisIdx(-1);
1002 for (auto i = 0; i < (G4int)pv->size(); ++i)
1003 {
1004 if ( (*pv)(i) == this ) { thisIdx = i; break; }
1005 }
1006 if ( thisIdx < 0 ) return false; // -- to ignore pure along processes
1007 for (std::size_t i=0; i<(fSharedData->fBiasingProcessInterfaces).size(); ++i)
1008 {
1009 if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
1010 {
1011 G4int thatIdx(-1);
1012 for (auto j = 0; j < (G4int)pv->size(); ++j)
1013 {
1014 if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] )
1015 {
1016 thatIdx = j; break;
1017 }
1018 }
1019 if ( thatIdx >= 0 ) // -- to ignore pure along processes
1020 {
1021 if ( thisIdx < thatIdx )
1022 {
1023 isLast = false;
1024 break;
1025 }
1026 }
1027 }
1028 }
1029 return isLast;
1030}
1031
1032void G4BiasingProcessInterface::SetUpFirstLastFlags()
1033{
1034 for (G4int iPhys = 0; iPhys < 2; ++iPhys)
1035 {
1036 G4bool physOnly = ( iPhys == 1 );
1037 fFirstLastFlags[IdxFirstLast( 1, 1, iPhys)] = IsFirstPostStepGPILInterface(physOnly);
1038 fFirstLastFlags[IdxFirstLast( 0, 1, iPhys)] = IsLastPostStepGPILInterface(physOnly);
1039 fFirstLastFlags[IdxFirstLast( 1, 0, iPhys)] = IsFirstPostStepDoItInterface(physOnly);
1040 fFirstLastFlags[IdxFirstLast( 0, 0, iPhys)] = IsLastPostStepDoItInterface(physOnly);
1041 }
1042
1043 // -- for itself, for optimization:
1044 fIamFirstGPIL = GetIsFirstPostStepGPILInterface( false );
1045}
1046
1047void G4BiasingProcessInterface::ResetForUnbiasedTracking()
1048{
1049 fOccurenceBiasingOperation = nullptr;
1050 fFinalStateBiasingOperation = nullptr;
1051 fNonPhysicsBiasingOperation = nullptr;
1052 fBiasingInteractionLaw = nullptr;
1053}
1054
1055void G4BiasingProcessInterface::
1056InvokeWrappedProcessPostStepGPIL( const G4Track& track,
1057 G4double previousStepSize,
1059{
1060 G4double usedPreviousStepSize = previousStepSize;
1061 // -- if the physics process has been under occurrence biasing in the previous step
1062 // -- we reset it, as we don't know if it will be biased again or not in this
1063 // -- step. The pity is that PostStepGPIL and interaction length (cross-section)
1064 // -- calculations are done both in the PostStepGPIL of the process, while here we
1065 // -- are just interested in the calculation of the cross-section. This is a pity
1066 // -- as this forces to re-generated a random number for nothing.
1067 if ( fResetWrappedProcessInteractionLength )
1068 {
1069 fResetWrappedProcessInteractionLength = false;
1070 fWrappedProcess->ResetNumberOfInteractionLengthLeft();
1071 // -- We set "previous step size" as 0.0, to let the process believe this is first step:
1072 usedPreviousStepSize = 0.0;
1073 }
1074 // -- GPIL response:
1075 fWrappedProcessPostStepGPIL = fWrappedProcess->PostStepGetPhysicalInteractionLength(track, usedPreviousStepSize, condition);
1076 fWrappedProcessForceCondition = *condition;
1077 // -- and (inverse) cross-section:
1078 fWrappedProcessInteractionLength = fWrappedProcess->GetCurrentInteractionLength();
1079}
1080
1081void G4BiasingProcessInterface::ReorderBiasingVectorAsGPIL()
1082{
1083 // -- re-order vector of processes to match that of the GPIL:
1084 std::vector < G4BiasingProcessInterface* > tmpProcess ( fSharedData->fBiasingProcessInterfaces );
1085 ( fSharedData -> fBiasingProcessInterfaces ) . clear();
1086 ( fSharedData -> fPhysicsBiasingProcessInterfaces ) . clear();
1087 ( fSharedData -> fNonPhysicsBiasingProcessInterfaces ) . clear();
1088 ( fSharedData -> fPublicBiasingProcessInterfaces ) . clear();
1089 ( fSharedData -> fPublicPhysicsBiasingProcessInterfaces ) . clear();
1090 ( fSharedData -> fPublicNonPhysicsBiasingProcessInterfaces ) . clear();
1091
1092 const G4ProcessVector* pv = fProcessManager->GetPostStepProcessVector(typeGPIL);
1093 for (auto i = 0; i < (G4int)pv->size(); ++i)
1094 {
1095 for (std::size_t j = 0; j < tmpProcess.size(); ++j)
1096 {
1097 if ( (*pv)(i) == tmpProcess[j] )
1098 {
1099 ( fSharedData->fBiasingProcessInterfaces ).push_back( tmpProcess[j] );
1100 ( fSharedData->fPublicBiasingProcessInterfaces ).push_back( tmpProcess[j] );
1101 if ( tmpProcess[j] -> fIsPhysicsBasedBiasing )
1102 {
1103 ( fSharedData->fPhysicsBiasingProcessInterfaces ).push_back( tmpProcess[j] );
1104 ( fSharedData->fPublicPhysicsBiasingProcessInterfaces ).push_back( tmpProcess[j] );
1105 }
1106 else
1107 {
1108 ( fSharedData -> fNonPhysicsBiasingProcessInterfaces ).push_back( tmpProcess[j] );
1109 ( fSharedData -> fPublicNonPhysicsBiasingProcessInterfaces ).push_back( tmpProcess[j] );
1110 }
1111 break;
1112 }
1113 }
1114 }
1115}
G4BiasingAppliedCase
@ BAC_FinalState
@ BAC_NonPhysics
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ForceCondition
@ NotForced
G4GPILSelection
@ NotCandidateForSelection
@ typeGPIL
@ typeDoIt
@ fGeomBoundary
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4bool GetIsFirstPostStepGPILInterface(G4bool physOnly=true) const
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
G4bool IsLastPostStepGPILInterface(G4bool physOnly=true) const
virtual void SetProcessManager(const G4ProcessManager *)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual void SetMasterProcess(G4VProcess *masterP)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
G4bool GetIsLastPostStepDoItInterface(G4bool physOnly=true) const
virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition &pd)
G4bool GetIsLastPostStepGPILInterface(G4bool physOnly=true) const
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly=true) const
G4bool IsFirstPostStepDoItInterface(G4bool physOnly=true) const
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4bool IsLastPostStepDoItInterface(G4bool physOnly=true) const
G4BiasingProcessInterface(const G4String &name="biasWrapper(0)")
const G4BiasingProcessSharedData * GetSharedData() const
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *pd, const G4String &s, G4bool f)
virtual const G4ProcessManager * GetProcessManager()
virtual void PreparePhysicsTable(const G4ParticleDefinition &pd)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
virtual void BuildPhysicsTable(const G4ParticleDefinition &pd)
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *pd, const G4String &s, G4bool f)
virtual G4bool IsApplicable(const G4ParticleDefinition &pd)
G4bool IsFirstPostStepGPILInterface(G4bool physOnly=true) const
virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition &pd)
std::size_t size() const
G4double GetStepLength() const
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetVolume() const
G4int GetCurrentStepNumber() const
const G4Step * GetStep() const
const G4String & GetName() const
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators()
static G4VBiasingOperator * GetBiasingOperator(const G4LogicalVolume *)
G4TrackStatus GetTrackStatus() const
G4LogicalVolume * GetLogicalVolume() const
const G4VProcess * GetMasterProcess() const
virtual const G4ProcessManager * GetProcessManager()
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition G4VProcess.cc:46
G4ProcessType GetProcessType() const
virtual void SetMasterProcess(G4VProcess *masterP)
void SetProcessSubType(G4int)
virtual void SetProcessManager(const G4ProcessManager *)
const G4String & GetProcessName() const
#define DBL_MAX
Definition templates.hh:62