Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
GB05BOptrSplitAndKillByCrossSection.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 // $Id$
27 //
30 
34 
35 #include "G4ParticleDefinition.hh"
36 #include "G4ParticleTable.hh"
37 #include "G4VProcess.hh"
38 
39 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
40 
43  G4String name)
44  : G4VBiasingOperator(name),
45  fSetup(true)
46 {
47  fParticleToBias = G4ParticleTable::GetParticleTable()->FindParticle(particleName);
48 
49  if ( fParticleToBias == 0 )
50  {
52  ed << "Particle `" << particleName << "' not found !" << G4endl;
53  G4Exception("GB05BOptrSplitAndKillByCrossSection(...)",
54  "exGB05.01",
56  ed);
57  }
58 
59  fSplitAndKillByCrossSection =
60  new GB05BOptnSplitAndKillByCrossSection("splitterFor_"+particleName);
61 
62 }
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 
67 {
68  delete fSplitAndKillByCrossSection;
69 }
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
74 {
75  // ---------------
76  // -- Setup stage:
77  // ---------------
78  // -- Start by collecting the pointer of the physics processes
79  // -- considered for the splitting by cross-sections. Doing so,
80  // -- this also verifies that these physics processes are each
81  // -- under control of a G4BiasingProcessInterface wrapper.
82  if ( fSetup )
83  {
84  const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
85  const G4BiasingProcessSharedData* sharedData =
87  if ( sharedData )
88  {
89  for ( size_t i = 0 ; i < fProcessesToEquipoise.size() ; i++ )
90  {
91  G4bool processFound(false);
92  for ( size_t j = 0 ;
93  j < (sharedData->GetPhysicsBiasingProcessInterfaces()).size();
94  j++ )
95  {
96  const G4BiasingProcessInterface* wrapperProcess =
97  (sharedData->GetPhysicsBiasingProcessInterfaces())[j];
98  if ( fProcessesToEquipoise[i] ==
99  wrapperProcess->GetWrappedProcess()->GetProcessName() )
100  {
101  fProcesses.push_back( wrapperProcess->GetWrappedProcess() );
102  processFound = true;
103  break;
104  }
105  }
106  if ( !processFound )
107  {
108  G4String particleName = "(unknown)";
109  if ( fParticleToBias != nullptr )
110  {
111  particleName = fParticleToBias->GetParticleName();
112  }
114  ed << "Process `" << fProcessesToEquipoise[i]
115  << "' not found for particle `" << particleName << "'"
116  << G4endl;
117  G4Exception("GB05BOptrSplitAndKillByCrossSection::StartRun(...)",
118  "exGB05.02",
119  JustWarning,
120  ed);
121  }
122  }
123  }
124  fSetup = false;
125  }
126 
127  if ( fProcessesToEquipoise.size() == 0 || fProcesses.size() == 0 )
128  {
130  ed << "No processes to counterbalance for defined or found ! "
131  << "Biasing will do nothing."
132  << G4endl;
133  G4Exception("GB05BOptrSplitAndKillByCrossSection::StartRun(...)",
134  "exGB05.03",
135  JustWarning,
136  ed);
137  }
138 
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142 
144 GB05BOptrSplitAndKillByCrossSection::
145 ProposeNonPhysicsBiasingOperation(const G4Track* track,
146  const G4BiasingProcessInterface* /*callingProcess*/)
147 {
148 
149  // -----------------------------------------------------
150  // -- Check if current particle type is the one to bias:
151  // -----------------------------------------------------
152  if ( track->GetDefinition() != fParticleToBias ) return 0;
153 
154  // --------------------------------------------------------------------
155  // -- Compute the total cross-section for the physics processes
156  // -- considered.
157  // -- These physics processes have been updated to the current track
158  // -- state by their related wrapper G4BiasingProcessInterface objects,
159  // -- the cross-sections/mean free pathes are hence usable.
160  // --------------------------------------------------------------------
161  G4double totalCrossSection(0.0);
162  for ( size_t i = 0 ; i < fProcesses.size() ; i++ )
163  {
164  G4double interactionLength = fProcesses[i]->GetCurrentInteractionLength();
165  if ( interactionLength < DBL_MAX/10. )
166  totalCrossSection += 1./interactionLength;
167  }
168  if ( totalCrossSection < DBL_MIN ) return 0;
169 
170  G4double totalInteractionLength = 1./totalCrossSection;
171 
172  // ---------------------------------------------------------------------
173  // -- Passes the updated "absorption" cross-section (interaction length)
174  // -- to the biasing operation, and returns this operation:
175  // ---------------------------------------------------------------------
176  fSplitAndKillByCrossSection->SetInteractionLength( totalInteractionLength );
177 
178  return fSplitAndKillByCrossSection;
179 
180 }
181 
184 {
185  fProcessesToEquipoise.push_back( processName );
186 }
G4ParticleDefinition * GetDefinition() const
const XML_Char * name
Definition: expat.h:151
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const std::vector< const G4BiasingProcessInterface * > & GetPhysicsBiasingProcessInterfaces() const
Definition of the GB05BOptnSplitAndKillByCrossSection class.
const G4String & GetParticleName() const
G4VProcess * GetWrappedProcess() const
GB05BOptrSplitAndKillByCrossSection(G4String particleToBias, G4String name="SplitAndKillByXS")
bool G4bool
Definition: G4Types.hh:79
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4ParticleTable * GetParticleTable()
G4ProcessManager * GetProcessManager() const
#define DBL_MIN
Definition: templates.hh:75
Definition of the GB05BOptrSplitAndKillByCrossSection class.
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
const G4BiasingProcessSharedData * GetSharedData() const