Geant4  10.01.p02
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 //
27 #include "G4VBiasingOperator.hh"
28 #include "G4VBiasingOperation.hh"
30 #include "G4ParticleChange.hh"
34 #include "G4ProcessManager.hh"
35 #include "G4BiasingTrackData.hh"
37 #include "G4BiasingAppliedCase.hh"
38 
43 
44 
46  : G4VProcess ( name ),
47  fWrappedProcess ( 0 ),
48  fIsPhysicsBasedBiasing ( false ),
49  fWrappedProcessIsAtRest ( false ),
50  fWrappedProcessIsAlong ( false ),
51  fWrappedProcessIsPost ( false ),
52  fWrappedProcessInteractionLength( -1.0 ),
53  fBiasingInteractionLaw ( 0 ),
54  fPhysicalInteractionLaw ( 0 ),
55  fOccurenceBiasingParticleChange ( 0 ),
56  fIamFirstGPIL ( false ),
57  fSharedData ( 0 )
58 
59 {
60  for (G4int i = 0 ; i < 8 ; i++) fFirstLastFlags[i] = false;
61  fResetInteractionLaws.Put( true );
62  fCommonStart.Put(true);
63  fCommonEnd.Put(true);
64 }
65 
66 
68  G4bool wrappedIsAtRest, G4bool wrappedIsAlongStep, G4bool wrappedIsPostStep,
69  G4String useThisName)
70  : G4VProcess( useThisName != "" ? useThisName : "biasWrapper("+wrappedProcess->GetProcessName()+")",
71  wrappedProcess->GetProcessType()),
72  fWrappedProcess ( wrappedProcess ),
73  fIsPhysicsBasedBiasing ( true ),
74  fWrappedProcessIsAtRest ( wrappedIsAtRest ),
75  fWrappedProcessIsAlong ( wrappedIsAlongStep ),
76  fWrappedProcessIsPost ( wrappedIsPostStep ),
77  fWrappedProcessInteractionLength( -1.0 ),
78  fBiasingInteractionLaw ( 0 ),
79  fPhysicalInteractionLaw ( 0 ),
80  fOccurenceBiasingParticleChange ( 0 ),
81  fIamFirstGPIL ( false ),
82  fSharedData ( 0 )
83 {
84  for (G4int i = 0 ; i < 8 ; i++) fFirstLastFlags[i] = false;
85  fResetInteractionLaws.Put( true );
86  fCommonStart.Put(true);
87  fCommonEnd.Put(true);
88 
90 
91  // -- create physical interaction law:
92  fPhysicalInteractionLaw = new G4InteractionLawPhysical("PhysicalInteractionLawFor("+GetProcessName()+")");
93  // -- instantiate particle change wrapper for occurence biaising:
96  // -- instantiate a "do nothing" particle change:
98 }
99 
100 
101 
103 {
106  if ( fParticleChange ) delete fParticleChange;
108 }
109 
110 
112 {
113  G4MapCache< const G4ProcessManager*,
114  G4BiasingProcessSharedData* >::const_iterator itr = fSharedDataMap.Find( mgr );
115  if ( itr != fSharedDataMap.End( ) )
116  {
117  return (*itr).second;
118  }
119  else return 0;
120 }
121 
122 
124 {
125  fCurrentTrack = track;
135 
136  fPreviousStepSize = -1.0;
137 
139 
140  if ( fCommonStart.Get() )
141  {
142  fCommonStart.Put( false );// = false;
143  fCommonEnd.Put(true);// = true;
144 
145  fSharedData-> fCurrentBiasingOperator = 0;
147 
148  // -- §§ Add a "fSharedData->nStarting" here and outside bracket "fSharedData->nStarting++" and " if (fSharedData->nStarting) == fSharedData->(vector interface length)"
149  // -- §§ call to the loop "StartTracking" of operators"
150 
151  for ( size_t optr = 0 ; optr < ( G4VBiasingOperator::GetBiasingOperators() ).size() ; optr ++)
153  }
154 }
155 
156 
158 {
162 
163 
164  // -- !! this part might have to be improved : could be time consuming
165  // -- !! and assumes all tracks are killed during tracking, which is
166  // -- !! not true : stacking operations may kill tracks
167  if ( fCommonEnd.Get() )
168  {
169  fCommonEnd.Put( false );// = false;
170  fCommonStart.Put( true );// = true;
171 
172  for ( size_t optr = 0 ; optr < ( G4VBiasingOperator::GetBiasingOperators() ).size() ; optr ++)
174 
175  // -- §§ for above loop, do as in StartTracking.
176 
178  {
180  if ( biasingData ) delete biasingData; // -- this also deregisters the biasing data from the track data store
182  {
183  const G4TrackVector* secondaries = fCurrentTrack->GetStep()->GetSecondary();
184  for ( size_t i2nd = 0 ; i2nd < secondaries->size() ; i2nd++ )
185  {
186  biasingData = G4BiasingTrackDataStore::GetInstance()->GetBiasingTrackData( (*secondaries)[i2nd] );
187  if ( biasingData ) delete biasingData;
188  }
189  }
190  }
191  }
192 }
193 
194 
195 
197  G4double previousStepSize,
199 {
200  // ---------------------------------------------------------------------------------------------
201  // -- The "biasing master" takes care for all biasing processes of update of biasing operators
202  // -- and invokes all PostStepGPIL of physical wrapped processes (anticipate stepping manager
203  // -- call ! ) so that all cross-sections are updated with current step, and available right
204  // -- away to the biasing operator.
205  // ---------------------------------------------------------------------------------------------
206  if ( fIamFirstGPIL )
207  {
208  // -- Update previous biasing operator:
210  // -- If new volume, get possible new biasing operator:
211  // ----------------------------------------------------
212  // -- [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 ?]
213  G4bool firstStepInVolume = ( (track.GetStep()->GetPreStepPoint()->GetStepStatus() == fGeomBoundary) || (track.GetCurrentStepNumber() == 1) );
214  fSharedData->fIsNewOperator = false;
216  if ( firstStepInVolume )
217  {
219  fSharedData->fCurrentBiasingOperator = newOperator ;
220  if ( newOperator != fSharedData->fPreviousBiasingOperator )
221  {
223  fSharedData->fIsNewOperator = ( newOperator != 0 );
224  }
225  }
226 
227 
228  // -- calls to wrapped process PostStepGPIL's:
229  // -------------------------------------------
230  // -- Each physics wrapper process has its
231  // -- fWrappedProcessPostStepGPIL ,
232  // -- fWrappedProcessForceCondition ,
233  // -- fWrappedProcessInteractionLength
234  // -- updated.
236  {
237  for ( size_t i = 0 ; i < (fSharedData->fPhysicsBiasingProcessInterfaces).size(); i++ )
238  (fSharedData->fPhysicsBiasingProcessInterfaces)[i]->InvokeWrappedProcessPostStepGPIL( track, previousStepSize, condition );
239  }
240  }
241 
242 
243  // -- Remember previous operator and proposed operations, if any, and reset:
244  // -------------------------------------------------------------------------
245  // -- remember only in case some biasing might be called
246  if ( ( fSharedData->fPreviousBiasingOperator != 0 ) ||
248  {
253  // -- reset:
258  // -- Physics PostStep and AlongStep GPIL
259  // fWrappedProcessPostStepGPIL : updated by InvokeWrappedProcessPostStepGPIL(...) above
261  // fWrappedProcessInteractionLength : updated by InvokeWrappedProcessPostStepGPIL(...) above; inverse of analog cross-section.
262  // fWrappedProcessForceCondition : updated by InvokeWrappedProcessPostStepGPIL(...) above
268  // -- for helper:
269  fPreviousStepSize = previousStepSize;
270  }
271 
272 
273  // -- previous step size value; it is switched to zero if resetting a wrapped process:
274  // -- (same trick used than in InvokedWrappedProcessPostStepGPIL )
275  G4double usedPreviousStepSize = previousStepSize;
276 
277  // ----------------------------------------------
278  // -- If leaving a biasing operator, let it know:
279  // ----------------------------------------------
281  {
282  (fSharedData->fPreviousBiasingOperator)->ExitingBiasing( &track, this );
283  // -- if no further biasing operator, reset process behavior to standard tracking:
285  {
288  {
289  // -- if the physics process has been under occurence biasing, reset it:
291  {
294  // -- We set "previous step size" as 0.0, to let the process believe this is first step:
295  usedPreviousStepSize = 0.0;
296  }
297  }
298  }
299  }
300 
301 
302  // --------------------------------------------------------------
303  // -- no operator : analog tracking if physics-based, or nothing:
304  // --------------------------------------------------------------
306  {
307  // -- take note of the "usedPreviousStepSize" value:
308  if ( fIsPhysicsBasedBiasing ) return fWrappedProcess->PostStepGetPhysicalInteractionLength(track, usedPreviousStepSize, condition);
309  else
310  {
311  *condition = NotForced;
312  return DBL_MAX;
313  }
314  }
315 
316 
317 
318  // --------------------------------------------------
319  // -- An biasing operator exists. Proceed with
320  // -- treating non-physics and physics biasing cases:
321  //---------------------------------------------------
322 
323  // -- non-physics-based biasing case:
324  // ----------------------------------
325  if ( !fIsPhysicsBasedBiasing )
326  {
327  fNonPhysicsBiasingOperation = (fSharedData->fCurrentBiasingOperator)->GetProposedNonPhysicsBiasingOperation( &track, this );
328  if ( fNonPhysicsBiasingOperation == 0 )
329  {
330  *condition = NotForced;
331  return DBL_MAX;
332  }
333  return fNonPhysicsBiasingOperation->DistanceToApplyOperation(&track, previousStepSize, condition);
334  }
335 
336 
337  // -- Physics based biasing case:
338  // ------------------------------
339  // -- Ask for possible GPIL biasing operation:
340  fOccurenceBiasingOperation = (fSharedData->fCurrentBiasingOperator)->GetProposedOccurenceBiasingOperation( &track, this );
341 
342 
343  // -- no operation for occurence biasing, analog GPIL returns the wrapped process GPIL and condition values
344  if ( fOccurenceBiasingOperation == 0 )
345  {
346  *condition = fWrappedProcessForceCondition;
348  }
349 
350  // -- A valid GPIL biasing operation has been proposed:
351  // -- 0) remember wrapped process will need to be reset on biasing exit, if particle survives:
353  // -- 1) update process interaction length for reference analog interaction law ( fWrappedProcessInteractionLength updated/collected above):
355  // -- 2) Collect biasing interaction law:
356  // -- The interaction law pointer is collected as a const pointer to the interaction law object.
357  // -- This interaction law will be kept under control of the biasing operation, which is the only
358  // -- entity that will change the state of the biasing interaction law.
359  // -- The force condition for biasing is asked at the same time, passing the analog one as default:
362  // -- 3) Ask operation to sample the biasing interaction law:
364 
365  // -- finish
366  *condition = fBiasingForceCondition;
367  return fBiasingPostStepGPIL;
368 
369 }
370 
371 
372 
374  const G4Step& step)
375 {
376  // ---------------------------------------
377  // -- case outside of volume with biasing:
378  // ---------------------------------------
379  if ( fSharedData->fCurrentBiasingOperator == 0 ) return fWrappedProcess->PostStepDoIt(track, step);
380 
381  // ----------------------------
382  // -- non-physics biasing case:
383  // ----------------------------
384  if ( !fIsPhysicsBasedBiasing )
385  {
386  G4VParticleChange* particleChange = fNonPhysicsBiasingOperation->GenerateBiasingFinalState( &track, &step );
387  (fSharedData->fCurrentBiasingOperator)->ReportOperationApplied( this, BAC_NonPhysics, fNonPhysicsBiasingOperation, particleChange );
388  return particleChange;
389  }
390 
391  // -- physics biasing case:
392  // ------------------------
393  // -- It proceeds with the following logic:
394  // -- 1) Obtain the final state
395  // -- This final state may be analog or biased.
396  // -- The biased final state is obtained through a biasing operator
397  // -- returned by the operator.
398  // -- 2) The biased final state may be asked to be "force as it is"
399  // -- in what case the particle change is returned as is to the
400  // -- stepping.
401  // -- In all other cases (analog final state or biased final but
402  // -- not forced) the final state weight may be modified by the
403  // -- occurence biasing, if such an occurence biasing is at play.
404 
405  // -- Get final state, biased or analog:
406  G4VParticleChange* finalStateParticleChange;
408  fFinalStateBiasingOperation = (fSharedData->fCurrentBiasingOperator)->GetProposedFinalStateBiasingOperation( &track, this );
409  // -- Flag below is to force the biased generated particle change to be retruned "as is" to the stepping, disregarding there
410  // -- was or not a occurence biasing that would apply. Weight relevance under full responsibility of the biasing operation.
411  G4bool forceBiasedFinalState = false;
412  if ( fFinalStateBiasingOperation != 0 )
413  {
414  finalStateParticleChange = fFinalStateBiasingOperation->ApplyFinalStateBiasing( this, &track, &step, forceBiasedFinalState );
415  BAC = BAC_FinalState;
416  }
417  else
418  {
419  finalStateParticleChange = fWrappedProcess->PostStepDoIt(track, step);
420  BAC = BAC_None ;
421  }
422 
423  // -- if no occurence biasing operation, we're done:
424  if ( fOccurenceBiasingOperation == 0 )
425  {
426  (fSharedData->fCurrentBiasingOperator)->ReportOperationApplied( this, BAC, fFinalStateBiasingOperation, finalStateParticleChange );
427  return finalStateParticleChange;
428  }
429 
430  // -- if biased final state has been asked to be forced, we're done:
431  if ( forceBiasedFinalState )
432  {
433  (fSharedData->fCurrentBiasingOperator)->ReportOperationApplied( this, BAC, fFinalStateBiasingOperation, finalStateParticleChange );
434  return finalStateParticleChange;
435  }
436 
437  // -- If occurence biasing, applies the occurence biasing weight correction on top of final state (biased or not):
438  G4double weightForInteraction = 1.0;
439  if ( !fBiasingInteractionLaw->IsSingular() ) weightForInteraction =
441  fBiasingInteractionLaw ->ComputeEffectiveCrossSectionAt(step.GetStepLength());
442  else
443  {
444  // -- at this point effective XS can only be infinite, if not, there is a logic problem
446  {
448  ed << "Internal inconsistency in cross-section handling. Please report !" << G4endl;
449  G4Exception(" G4BiasingProcessInterface::PostStepDoIt(...)",
450  "BIAS.GEN.02",
451  JustWarning,
452  ed);
453  // -- if XS is infinite, weight is zero (and will stay zero), but we'll do differently.
454  // -- Should foresee in addition something to remember that in case of singular
455  // -- distribution, weight can only be partly calculated
456  }
457  }
458 
459  if ( weightForInteraction <= 0. )
460  {
462  ed << " Negative interaction weight : w_I = "
463  << weightForInteraction <<
466  " step length = " << step.GetStepLength() <<
467  " Interaction law = `" << fBiasingInteractionLaw << "'" <<
468  G4endl;
469  G4Exception(" G4BiasingProcessInterface::PostStepDoIt(...)",
470  "BIAS.GEN.03",
471  JustWarning,
472  ed);
473 
474  }
475 
476  (fSharedData->fCurrentBiasingOperator)->ReportOperationApplied( this, BAC,
477  fOccurenceBiasingOperation, weightForInteraction,
478  fFinalStateBiasingOperation, finalStateParticleChange );
479 
482  fOccurenceBiasingParticleChange->SetWrappedParticleChange( finalStateParticleChange );
483  fOccurenceBiasingParticleChange->ProposeTrackStatus( finalStateParticleChange->GetTrackStatus() );
484  fOccurenceBiasingParticleChange->StealSecondaries(); // -- this also makes weightForInteraction applied to secondaries stolen
485 
486  // -- finish:
488 
489 }
490 
491 
492 // -- AlongStep methods:
494  G4double previousStepSize,
495  G4double currentMinimumStep,
496  G4double& proposedSafety,
497  G4GPILSelection* selection)
498 {
499  // -- for helper methods:
500  fCurrentMinimumStep = currentMinimumStep;
501  fProposedSafety = proposedSafety;
502 
503 
504  // -- initialization default case:
506  *selection = NotCandidateForSelection;
507  // ---------------------------------------
508  // -- case outside of volume with biasing:
509  // ---------------------------------------
511  {
514  previousStepSize,
515  currentMinimumStep,
516  proposedSafety,
517  selection);
519  }
520 
521  // --------------------------------------------------------------------
522  // -- non-physics based biasing: no along operation expected (for now):
523  // --------------------------------------------------------------------
525 
526  // ----------------------
527  // -- physics-based case:
528  // ----------------------
529  if ( fOccurenceBiasingOperation == 0 )
530  {
533  previousStepSize,
534  currentMinimumStep,
535  proposedSafety,
536  selection);
538  }
539 
540 
541  // ----------------------------------------------------------
542  // -- From here we have an valid occurence biasing operation:
543  // ----------------------------------------------------------
544  // -- Give operation opportunity to shorten step proposed by physics process:
546  G4double minimumStep = fBiasingAlongStepGPIL < currentMinimumStep ? fBiasingAlongStepGPIL : currentMinimumStep ;
547  // -- wrapped process is called with minimum step ( <= currentMinimumStep passed ) : an along process can not
548  // -- have its operation stretched over what it expects:
550  {
552  previousStepSize,
553  minimumStep,
554  proposedSafety,
555  selection);
556  fWrappedProcessGPILSelection = *selection;
558  }
559  else
560  {
563  }
564 
565  *selection = fBiasingGPILSelection;
566 
568 
569 }
570 
572  const G4Step& step)
573 {
574  // ---------------------------------------
575  // -- case outside of volume with biasing:
576  // ---------------------------------------
578  {
579  if ( fWrappedProcessIsAlong ) return fWrappedProcess->AlongStepDoIt(track, step);
580  else
581  {
583  return fDummyParticleChange;
584  }
585  }
586 
587  // -----------------------------------
588  // -- case inside volume with biasing:
589  // -----------------------------------
591  else
592  {
595  }
596  G4double weightForNonInteraction (1.0);
597  if ( fBiasingInteractionLaw != 0 )
598  {
599  weightForNonInteraction =
601  fBiasingInteractionLaw ->ComputeNonInteractionProbabilityAt(step.GetStepLength());
602 
603  fOccurenceBiasingOperation->AlongMoveBy( this, &step, weightForNonInteraction );
604 
605  if ( weightForNonInteraction <= 0. )
606  {
608  ed << " Negative non interaction weight : w_NI = " << weightForNonInteraction <<
610  " p_NI(bias) = " << fBiasingInteractionLaw ->ComputeNonInteractionProbabilityAt(step.GetStepLength()) <<
611  " step length = " << step.GetStepLength() <<
612  " biasing interaction law = `" << fBiasingInteractionLaw->GetName() << "'" << G4endl;
613  G4Exception(" G4BiasingProcessInterface::AlongStepDoIt(...)",
614  "BIAS.GEN.04",
615  JustWarning,
616  ed);
617  }
618 
619  }
620 
622 
624 
625 }
626 
627 // -- AtRest methods
630 {
631  return fWrappedProcess->AtRestGetPhysicalInteractionLength(track, condition);
632 }
634  const G4Step& step)
635 {
636  return fWrappedProcess->AtRestDoIt(track, step);
637 }
638 
639 
641 {
642  if ( fWrappedProcess != 0 ) return fWrappedProcess->IsApplicable(pd);
643  else return true;
644 }
645 
646 
648 {
649  // -- Master for this process:
651  // -- Master for wrapped process:
652  if ( fWrappedProcess != 0 )
653  {
654  const G4BiasingProcessInterface* thisWrapperMaster = (const G4BiasingProcessInterface *)GetMasterProcess();
655  // -- paranoia check: (?)
656  G4VProcess* wrappedMaster = 0;
657  wrappedMaster = thisWrapperMaster->GetWrappedProcess();
658  fWrappedProcess->SetMasterProcess( wrappedMaster );
659  }
660 }
661 
662 
664 {
665  // -- PreparePhysicsTable(...) has been called first for all processes,
666  // -- so the first/last flags and G4BiasingProcessInterface vector of processes have
667  // -- been properly setup, fIamFirstGPIL is valid.
668  if ( fWrappedProcess != 0 )
669  {
671  }
672 
673  // -- Re-order vector of processes to match that of the GPIL
674  // -- (made for fIamFirstGPIL, but important is to have it made once):
676 
677 }
678 
679 
681 {
682  // -- Let process finding its first/last position in the process manager:
684  if ( fWrappedProcess != 0 )
685  {
687  }
688 }
689 
690 
692 {
693  if ( fWrappedProcess != 0 ) return fWrappedProcess->StorePhysicsTable(pd, s, f);
694  else return false;
695 }
696 
697 
699 {
700  if ( fWrappedProcess != 0 ) return fWrappedProcess->RetrievePhysicsTable(pd, s, f);
701  else return false;
702 }
703 
704 
706 {
709 
710  // -- initialize fSharedData pointer:
711  if ( fSharedDataMap.Find(mgr) == fSharedDataMap.End() )
712  {
715  }
716  else fSharedData = fSharedDataMap[mgr] ;
717  // -- augment list of co-operating processes:
718  fSharedData-> fBiasingProcessInterfaces.push_back( this );
719  fSharedData-> fPublicBiasingProcessInterfaces.push_back( this );
720  if ( fIsPhysicsBasedBiasing )
721  {
722  fSharedData-> fPhysicsBiasingProcessInterfaces.push_back( this );
723  fSharedData-> fPublicPhysicsBiasingProcessInterfaces.push_back( this );
724  }
725  else
726  {
727  fSharedData-> fNonPhysicsBiasingProcessInterfaces.push_back( this );
728  fSharedData-> fPublicNonPhysicsBiasingProcessInterfaces.push_back( this );
729  }
730  // -- remember process manager:
731  fProcessManager = mgr;
732 }
733 
734 
736 {
737  if ( fWrappedProcess != 0 ) return fWrappedProcess->GetProcessManager();
738  else return G4VProcess::GetProcessManager();
739 }
740 
741 
743 {
744  // -- PrepareWorkerPhysicsTable(...) has been called first for all processes,
745  // -- so the first/last flags and G4BiasingProcessInterface vector of processes have
746  // -- been properly setup, fIamFirstGPIL is valid.
747  if ( fWrappedProcess != 0 )
748  {
750  }
751 
752  // -- Re-order vector of processes to match that of the GPIL
753  // -- (made for fIamFirstGPIL, but important is to have it made once):
755 }
756 
757 
759 {
760  // -- Let process finding its first/last position in the process manager:
762 
763  if ( fWrappedProcess != 0 )
764  {
766  }
767 }
768 
769 
771 {
773 }
774 
775 
777 {
778  G4int iPhys = ( physOnly ) ? 1 : 0;
779  return fFirstLastFlags[IdxFirstLast( 1, 1, iPhys)];
780 }
781 
782 
784 {
785  G4int iPhys = ( physOnly ) ? 1 : 0;
786  return fFirstLastFlags[IdxFirstLast( 0, 1, iPhys)];
787 }
788 
789 
791 {
792  G4int iPhys = ( physOnly ) ? 1 : 0;
793  return fFirstLastFlags[IdxFirstLast( 1, 0, iPhys)];
794 }
795 
796 
798 {
799  G4int iPhys = ( physOnly ) ? 1 : 0;
800  return fFirstLastFlags[IdxFirstLast( 0, 0, iPhys)];
801 }
802 
803 
805 {
806  G4bool isFirst = true;
808  G4int thisIdx(-1);
809  for (G4int i = 0; i < pv->size(); i++ ) if ( (*pv)(i) == this ) { thisIdx = i; break; }
810  // for ( size_t i = 0; i < fCoInterfaces->size(); i++ )
811  for ( size_t i = 0; i < (fSharedData->fBiasingProcessInterfaces).size(); i++ )
812  {
813  // if ( ( (*fCoInterfaces)[i]->GetWrappedProcess() != 0 ) || !physOnly )
814  if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
815  {
816  G4int thatIdx(-1);
817  // for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (*fCoInterfaces)[i] ) { thatIdx = j; break; }
818  for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] ) { thatIdx = j; break; }
819  if ( thisIdx > thatIdx )
820  {
821  isFirst = false;
822  break;
823  }
824  }
825  }
826  return isFirst;
827 }
828 
829 
831 {
832  G4bool isLast = true;
834  G4int thisIdx(-1);
835  for (G4int i = 0; i < pv->size(); i++ ) if ( (*pv)(i) == this ) { thisIdx = i; break; }
836  // for ( size_t i = 0; i < fCoInterfaces->size(); i++ )
837  for ( size_t i = 0; i < (fSharedData->fBiasingProcessInterfaces).size(); i++ )
838  {
839  // if ( ( (*fCoInterfaces)[i]->GetWrappedProcess() != 0 ) || !physOnly )
840  if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
841  {
842  G4int thatIdx(-1);
843  // for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (*fCoInterfaces)[i] ) { thatIdx = j; break; }
844  for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] ) { thatIdx = j; break; }
845  if ( thisIdx < thatIdx )
846  {
847  isLast = false;
848  break;
849  }
850  }
851  }
852  return isLast;
853 }
854 
855 
857 {
858  G4bool isFirst = true;
860  G4int thisIdx(-1);
861  for (G4int i = 0; i < pv->size(); i++ ) if ( (*pv)(i) == this ) { thisIdx = i; break; }
862  // for ( size_t i = 0; i < fCoInterfaces->size(); i++ )
863  for ( size_t i = 0; i < (fSharedData->fBiasingProcessInterfaces).size(); i++ )
864  {
865  // if ( ( (*fCoInterfaces)[i]->GetWrappedProcess() != 0 ) || !physOnly )
866  if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
867  {
868  G4int thatIdx(-1);
869  // for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (*fCoInterfaces)[i] ) { thatIdx = j; break; }
870  for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] ) { thatIdx = j; break; }
871  if ( thisIdx > thatIdx )
872  {
873  isFirst = false;
874  break;
875  }
876  }
877  }
878  return isFirst;
879 }
880 
881 
883 {
884  G4bool isLast = true;
886  G4int thisIdx(-1);
887  for (G4int i = 0; i < pv->size(); i++ ) if ( (*pv)(i) == this ) { thisIdx = i; break; }
888  // for ( size_t i = 0; i < fCoInterfaces->size(); i++ )
889  for ( size_t i = 0; i < (fSharedData->fBiasingProcessInterfaces).size(); i++ )
890  {
891  // if ( ( (*fCoInterfaces)[i]->GetWrappedProcess() != 0 ) || !physOnly )
892  if ( (fSharedData->fBiasingProcessInterfaces)[i]->fIsPhysicsBasedBiasing || !physOnly )
893  {
894  G4int thatIdx(-1);
895  // for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (*fCoInterfaces)[i] ) { thatIdx = j; break; }
896  for (G4int j = 0; j < pv->size(); j++ ) if ( (*pv)(j) == (fSharedData->fBiasingProcessInterfaces)[i] ) { thatIdx = j; break; }
897  if ( thisIdx < thatIdx )
898  {
899  isLast = false;
900  break;
901  }
902  }
903  }
904  return isLast;
905 }
906 
907 
909 {
910  for ( G4int iPhys = 0; iPhys < 2; iPhys++ )
911  {
912  G4bool physOnly = ( iPhys == 1 );
913  fFirstLastFlags[IdxFirstLast( 1, 1, iPhys)] = IsFirstPostStepGPILInterface(physOnly);
914  fFirstLastFlags[IdxFirstLast( 0, 1, iPhys)] = IsLastPostStepGPILInterface(physOnly);
915  fFirstLastFlags[IdxFirstLast( 1, 0, iPhys)] = IsFirstPostStepDoItInterface(physOnly);
916  fFirstLastFlags[IdxFirstLast( 0, 0, iPhys)] = IsLastPostStepDoItInterface(physOnly);
917  }
918 
919  // -- for itself, for optimization:
921 }
922 
923 
925 {
930 }
931 
932 
934  G4double previousStepSize,
936 {
937  G4double usedPreviousStepSize = previousStepSize;
938  // -- if the physics process has been under occurence biasing in the previous step
939  // -- we reset it, as we don't know if it will be biased again or not in this
940  // -- step. The pity is that PostStepGPIL and interaction length (cross-section)
941  // -- calculations are done both in the PostStepGPIL of the process, while here we
942  // -- are just interested in the calculation of the cross-section. This is a pity
943  // -- as this forces to re-generated a random number for nothing.
945  {
948  // -- We set "previous step size" as 0.0, to let the process believe this is first step:
949  usedPreviousStepSize = 0.0;
950  }
951  // -- GPIL response:
952  fWrappedProcessPostStepGPIL = fWrappedProcess->PostStepGetPhysicalInteractionLength(track, usedPreviousStepSize, condition);
954  // -- and (inverse) cross-section:
956 }
957 
958 
960 {
961  // -- re-order vector of processes to match that of the GPIL:
962  std::vector < G4BiasingProcessInterface* > tmpProcess ( fSharedData->fBiasingProcessInterfaces );
963  ( fSharedData -> fBiasingProcessInterfaces ) . clear();
964  ( fSharedData -> fPhysicsBiasingProcessInterfaces ) . clear();
965  ( fSharedData -> fNonPhysicsBiasingProcessInterfaces ) . clear();
966  ( fSharedData -> fPublicBiasingProcessInterfaces ) . clear();
967  ( fSharedData -> fPublicPhysicsBiasingProcessInterfaces ) . clear();
968  ( fSharedData -> fPublicNonPhysicsBiasingProcessInterfaces ) . clear();
969 
971  for (G4int i = 0; i < pv->size(); i++ )
972  {
973  for ( size_t j = 0; j < tmpProcess.size(); j++ )
974  {
975  if ( (*pv)(i) == tmpProcess[j] )
976  {
977  ( fSharedData -> fBiasingProcessInterfaces ) . push_back( tmpProcess[j] );
978  ( fSharedData -> fPublicBiasingProcessInterfaces ) . push_back( tmpProcess[j] );
979  if ( tmpProcess[j] -> fIsPhysicsBasedBiasing )
980  {
981  ( fSharedData -> fPhysicsBiasingProcessInterfaces ) . push_back( tmpProcess[j] );
982  ( fSharedData -> fPublicPhysicsBiasingProcessInterfaces ) . push_back( tmpProcess[j] );
983  }
984  else
985  {
986  ( fSharedData -> fNonPhysicsBiasingProcessInterfaces ) . push_back( tmpProcess[j] );
987  ( fSharedData -> fPublicNonPhysicsBiasingProcessInterfaces ) . push_back( tmpProcess[j] );
988  }
989  break;
990  }
991  }
992  }
993 }
virtual void SetMasterProcess(G4VProcess *masterP)
Definition: G4VProcess.cc:212
virtual void Initialize(const G4Track &track)
G4double condition(const G4ErrorSymMatrix &m)
static G4MapCache< const G4ProcessManager *, G4BiasingProcessSharedData * > fSharedDataMap
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition &pd)
G4VBiasingOperation * fPreviousFinalStateBiasingOperation
G4bool IsLastPostStepGPILInterface(G4bool physOnly=true) const
virtual void PreparePhysicsTable(const G4ParticleDefinition &pd)
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
virtual G4double ComputeNonInteractionProbabilityAt(G4double length) const
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)=0
virtual void SetProcessManager(const G4ProcessManager *)
Definition: G4VProcess.hh:508
G4VBiasingOperation * fFinalStateBiasingOperation
static G4Cache< G4bool > fCommonEnd
G4VBiasingOperation * fNonPhysicsBiasingOperation
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
static G4BiasingTrackDataStore * GetInstance()
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *pd, const G4String &s, G4bool f)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
iterator End()
Definition: G4Cache.hh:411
static G4Cache< G4bool > fResetInteractionLaws
G4double GetStepLength() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &pd)
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)=0
virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition &pd)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
G4String name
Definition: TRTMaterials.hh:40
G4StepStatus GetStepStatus() const
G4ParticleChangeForNothing * fDummyParticleChange
value_type & Get() const
Definition: G4Cache.hh:282
void InvokeWrappedProcessPostStepGPIL(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4TrackStatus GetTrackStatus() const
std::vector< G4BiasingProcessInterface * > fBiasingProcessInterfaces
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
virtual void SetProcessManager(const G4ProcessManager *)
virtual G4double ComputeEffectiveCrossSectionAt(G4double length) const =0
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
G4VBiasingOperation * fPreviousOccurenceBiasingOperation
const G4ProcessManager * fProcessManager
const G4Step * GetStep() const
int G4int
Definition: G4Types.hh:78
G4BiasingTrackData * GetBiasingTrackData(const G4Track *track)
G4BiasingAppliedCase
const G4VBiasingInteractionLaw * fPreviousBiasingInteractionLaw
virtual G4double ProposeAlongStepLimit(const G4BiasingProcessInterface *)
G4VProcess * GetWrappedProcess() const
virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition &part)
Definition: G4VProcess.cc:202
G4bool GetIsFirstPostStepGPILInterface(G4bool physOnly=true) const
static const double s
Definition: G4SIunits.hh:150
G4StepPoint * GetPreStepPoint() const
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:101
void SetSecondaryWeightByProcess(G4bool)
static G4Cache< G4bool > fCommonStart
std::vector< G4BiasingProcessInterface * > fPhysicsBiasingProcessInterfaces
G4int GetCurrentStepNumber() const
G4bool GetIsLastPostStepDoItInterface(G4bool physOnly=true) const
const G4VBiasingInteractionLaw * fBiasingInteractionLaw
G4bool GetIsLastPostStepGPILInterface(G4bool physOnly=true) const
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:217
virtual G4double ComputeEffectiveCrossSectionAt(G4double length) const
bool G4bool
Definition: G4Types.hh:79
G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly=true) const
G4VBiasingOperation * fPreviousNonPhysicsBiasingOperation
G4double GetCurrentInteractionLength() const
Definition: G4VProcess.hh:462
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
void SetPhysicalCrossSection(G4double crossSection)
Definition: G4Step.hh:76
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:236
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *)=0
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4bool IsApplicable(const G4ParticleDefinition &pd)
virtual const G4ProcessManager * GetProcessManager()
static G4VBiasingOperator * GetBiasingOperator(const G4LogicalVolume *)
virtual G4bool IsSingular() const
G4int size() const
virtual const G4ProcessManager * GetProcessManager()
Definition: G4VProcess.hh:514
iterator Find(const key_type &k)
Definition: G4Cache.hh:417
G4LogicalVolume * GetLogicalVolume() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:210
G4VBiasingOperation * fOccurenceBiasingOperation
std::vector< G4Track * > G4TrackVector
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0
G4ParticleChangeForOccurenceBiasing * fOccurenceBiasingParticleChange
G4bool IsFirstPostStepGPILInterface(G4bool physOnly=true) const
virtual void SetMasterProcess(G4VProcess *masterP)
G4BiasingProcessInterface(G4String name="biasWrapper(0)")
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *pd, const G4String &s, G4bool f)
virtual G4GPILSelection ProposeGPILSelection(const G4GPILSelection wrappedProcessSelection)
G4VPhysicalVolume * GetVolume() const
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
G4double GetSampledInteractionLength() const
#define G4endl
Definition: G4ios.hh:61
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators()
const G4String & GetName() const
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:231
virtual G4bool IsEffectiveCrossSectionInfinite() const
G4TrackStatus GetTrackStatus() const
G4bool IsFirstPostStepDoItInterface(G4bool physOnly=true) const
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &)=0
virtual void AlongMoveBy(const G4BiasingProcessInterface *, const G4Step *, G4double)
G4InteractionLawPhysical * fPhysicalInteractionLaw
G4BiasingProcessSharedData * fSharedData
virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.cc:207
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
virtual void EndTracking()
Definition: G4VProcess.cc:113
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
G4ForceCondition
G4bool IsLastPostStepDoItInterface(G4bool physOnly=true) const
G4int IdxFirstLast(G4int firstLast, G4int GPILDoIt, G4int physAll) const
#define DBL_MAX
Definition: templates.hh:83
G4VBiasingOperator * fPreviousBiasingOperator
void Put(const value_type &val) const
Definition: G4Cache.hh:286
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *)=0
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:205
G4VBiasingOperator * fCurrentBiasingOperator
const G4TrackVector * GetSecondary() const
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
G4GPILSelection
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
const G4BiasingProcessSharedData * GetSharedData() const