Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
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 new G4ChannelingOptrChangeCrossSection(particleName);
64 fParticlesToBias.push_back( particle );
65 fBOptrForParticle[ particle ] = optr;
66}
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69
71 G4ParticleTable::G4PTblDicIterator* aParticleIterator =
72 (G4ParticleTable::GetParticleTable())->GetIterator();
73
74 aParticleIterator->reset();
75
76 while( (*aParticleIterator)() ){
77 G4ParticleDefinition* particle = aParticleIterator->value();
78 if (particle->GetPDGCharge() !=0) {
79 AddParticle(particle->GetParticleName());
80 }
81 }
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85
87G4ChannelingOptrMultiParticleChangeCrossSection::
88ProposeOccurenceBiasingOperation(const G4Track* track,
89 const G4BiasingProcessInterface* callingProcess){
90
91 if ( fCurrentOperator ) return fCurrentOperator->
92 GetProposedOccurenceBiasingOperation(track, callingProcess);
93 else return 0;
94}
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97
99 const G4ParticleDefinition* definition = track->GetParticleDefinition();
100 std::map < const G4ParticleDefinition*, G4ChannelingOptrChangeCrossSection* > :: iterator
101 it = fBOptrForParticle.find( definition );
102 fCurrentOperator = 0;
103 if ( it != fBOptrForParticle.end() ) fCurrentOperator = (*it).second;
104 fnInteractions = 0;
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
109void
110G4ChannelingOptrMultiParticleChangeCrossSection::
111OperationApplied( const G4BiasingProcessInterface* callingProcess,
112 G4BiasingAppliedCase biasingCase,
113 G4VBiasingOperation* occurenceOperationApplied,
114 G4double weightForOccurenceInteraction,
115 G4VBiasingOperation* finalStateOperationApplied,
116 const G4VParticleChange* particleChangeProduced ){
117 fnInteractions++;
118 if ( fCurrentOperator ) fCurrentOperator->ReportOperationApplied( callingProcess,
119 biasingCase,
120 occurenceOperationApplied,
121 weightForOccurenceInteraction,
122 finalStateOperationApplied,
123 particleChangeProduced );
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4BiasingAppliedCase
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
const G4String & GetParticleName() const
void reset(G4bool ifSkipIon=true)
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4ParticleTableIterator< G4String, G4ParticleDefinition * > G4PTblDicIterator
const G4ParticleDefinition * GetParticleDefinition() const
G4VBiasingOperation * GetProposedOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
G4VBiasingOperator(const G4String &name)
void ReportOperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)