Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ChannelingOptrMultiParticleChangeCrossSection.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//
29
30#include "G4ProcessManager.hh"
32#include "G4ParticleTable.hh"
33
34#include "G4SystemOfUnits.hh"
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
37
39G4VBiasingOperator("ChannelingChangeXS-Many"),
40fCurrentOperator(0),
41fnInteractions(0){
43}
44
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46
48 const G4ParticleDefinition* particle =
50
51 if ( particle == 0 )
52 {
54 ed << "Particle `" << particleName << "' not found !" << G4endl;
55 G4Exception("G4ChannelingOptrMultiParticleChangeCrossSection::AddParticle(...)",
56 "G4Channeling",
58 ed);
59 return;
60 }
61
63 fParticlesToBias.push_back( particle );
64 fBOptrForParticle[ particle ] = optr;
65}
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
70 G4ParticleTable::G4PTblDicIterator* aParticleIterator =
71 (G4ParticleTable::GetParticleTable())->GetIterator();
72
73 aParticleIterator->reset();
74
75 while( (*aParticleIterator)() ){
76 G4ParticleDefinition* particle = aParticleIterator->value();
77 if (particle->GetPDGCharge() !=0) {
78 AddParticle(particle->GetParticleName());
79 }
80 }
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
86G4ChannelingOptrMultiParticleChangeCrossSection::
87ProposeOccurenceBiasingOperation(const G4Track* track,
88 const G4BiasingProcessInterface* callingProcess){
89
90 if ( fCurrentOperator ) return fCurrentOperator->
91 GetProposedOccurenceBiasingOperation(track, callingProcess);
92 else return 0;
93}
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96
98 const G4ParticleDefinition* definition = track->GetParticleDefinition();
99 std::map < const G4ParticleDefinition*, G4ChannelingOptrChangeCrossSection* > :: iterator
100 it = fBOptrForParticle.find( definition );
101 fCurrentOperator = 0;
102 if ( it != fBOptrForParticle.end() ) fCurrentOperator = (*it).second;
103 fnInteractions = 0;
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107
108void
109G4ChannelingOptrMultiParticleChangeCrossSection::
110OperationApplied( const G4BiasingProcessInterface* callingProcess,
111 G4BiasingAppliedCase biasingCase,
112 G4VBiasingOperation* occurenceOperationApplied,
113 G4double weightForOccurenceInteraction,
114 G4VBiasingOperation* finalStateOperationApplied,
115 const G4VParticleChange* particleChangeProduced ){
116 fnInteractions++;
117 if ( fCurrentOperator ) fCurrentOperator->ReportOperationApplied( callingProcess,
118 biasingCase,
119 occurenceOperationApplied,
120 weightForOccurenceInteraction,
121 finalStateOperationApplied,
122 particleChangeProduced );
123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4BiasingAppliedCase
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4double GetPDGCharge() const
const G4String & GetParticleName() const
void reset(G4bool ifSkipIon=true)
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
const G4ParticleDefinition * GetParticleDefinition() const
G4VBiasingOperation * GetProposedOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
void ReportOperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)