Loading...
Searching...
No Matches
GB01BOptrMultiParticleChangeCrossSection.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/// \file GB01/src/GB01BOptrMultiParticleChangeCrossSection.cc
27/// \brief Implementation of the GB01BOptrMultiParticleChangeCrossSection class
28//
30#include "G4BiasingProcessInterface.hh"
31
33#include "G4ParticleDefinition.hh"
34#include "G4ParticleTable.hh"
35
36#include "G4SystemOfUnits.hh"
37
38//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
39
43
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45
47{
48 const G4ParticleDefinition* particle =
49 G4ParticleTable::GetParticleTable()->FindParticle( particleName );
50
51 if ( particle == 0 )
52 {
53 G4ExceptionDescription ed;
54 ed << "Particle `" << particleName << "' not found !" << G4endl;
55 G4Exception("GB01BOptrMultiParticleChangeCrossSection::AddParticle(...)",
56 "exGB01.02",
57 JustWarning,
58 ed);
59 return;
60 }
61
63 fParticlesToBias.push_back( particle );
64 fBOptrForParticle[ particle ] = optr;
65}
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
72 const G4BiasingProcessInterface* callingProcess)
73{
74 // -- examples of limitations imposed to apply the biasing:
75 // -- limit application of biasing to primary particles only:
76 if ( track->GetParentID() != 0 ) return 0;
77 // -- limit to at most 5 biased interactions:
78 if ( fnInteractions > 4 ) return 0;
79 // -- and limit to a weight of at least 0.05:
80 if ( track->GetWeight() < 0.05 ) return 0;
81
83 GetProposedOccurenceBiasingOperation(track, callingProcess);
84 else return 0;
85}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88
90{
91 // -- fetch the underneath biasing operator, if any, for the current particle type:
92 const G4ParticleDefinition* definition = track->GetParticleDefinition();
93 std::map < const G4ParticleDefinition*, GB01BOptrChangeCrossSection* > :: iterator
94 it = fBOptrForParticle.find( definition );
96 if ( it != fBOptrForParticle.end() ) fCurrentOperator = (*it).second;
97
98 // -- reset count for number of biased interactions:
100}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103
104void
106OperationApplied( const G4BiasingProcessInterface* callingProcess,
107 G4BiasingAppliedCase biasingCase,
108 G4VBiasingOperation* occurenceOperationApplied,
109 G4double weightForOccurenceInteraction,
110 G4VBiasingOperation* finalStateOperationApplied,
111 const G4VParticleChange* particleChangeProduced )
112{
113 // -- count number of biased interactions:
115
116 // -- inform the underneath biasing operator that a biased interaction occured:
117 if ( fCurrentOperator ) fCurrentOperator->ReportOperationApplied( callingProcess,
118 biasingCase,
119 occurenceOperationApplied,
120 weightForOccurenceInteraction,
121 finalStateOperationApplied,
122 particleChangeProduced );
123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Definition of the GB01BOptrChangeCrossSection class.
Definition of the GB01BOptrMultiParticleChangeCrossSection class.
std::vector< const G4ParticleDefinition * > fParticlesToBias
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
std::map< const G4ParticleDefinition *, GB01BOptrChangeCrossSection * > fBOptrForParticle

Applications | User Support | Publications | Collaboration