Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChannelingOptrChangeCrossSection.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
31#include "G4ParticleTable.hh"
32#include "G4VProcess.hh"
33
34#include "Randomize.hh"
35
37
39#include "G4EmProcessSubType.hh"
40
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42
44 G4String name)
46fChannelingID(-1),
47fSetup(true){
48 fParticleToBias = G4ParticleTable::GetParticleTable()->FindParticle(particleName);
49
50 if ( fParticleToBias == 0 )
51 {
53 ed << "Particle `" << particleName << "' not found !" << G4endl;
54 G4Exception("G4ChannelingOptrChangeCrossSection(...)",
55 "G4Channeling",
57 ed);
58 }
59
60 fProcessToDensity["channeling"] = fDensityRatioNone;
61}
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
66 for ( std::map< const G4BiasingProcessInterface*, G4BOptnChangeCrossSection* >::iterator
67 it = fChangeCrossSectionOperations.begin() ;
68 it != fChangeCrossSectionOperations.end() ;
69 it++ ) delete (*it).second;
70}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
75 if ( fSetup ){
76 const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
77 const G4BiasingProcessSharedData* sharedData =
79 if ( sharedData ){
80 for ( size_t i = 0 ; i < (sharedData->GetPhysicsBiasingProcessInterfaces()).size(); i++ ){
81 const G4BiasingProcessInterface* wrapperProcess =
82 (sharedData->GetPhysicsBiasingProcessInterfaces())[i];
83 G4String processName = wrapperProcess->GetWrappedProcess()->GetProcessName();
84 G4String operationName = "channelingChangeXS-" + processName;
85 fChangeCrossSectionOperations[wrapperProcess] =
86 new G4BOptnChangeCrossSection(operationName);
87
88 G4ProcessType type = wrapperProcess->GetWrappedProcess()->GetProcessType();
89 G4int subType = wrapperProcess->GetWrappedProcess()->GetProcessSubType();
90
91 switch (type) {
92 case fNotDefined:
93 fProcessToDensity[processName] = fDensityRatioNotDefined;
94 break;
95 case fTransportation:
96 fProcessToDensity[processName] = fDensityRatioNone;
97 break;
99 if(subType == fCoulombScattering ||
100 subType == fMultipleScattering){
101 fProcessToDensity[processName] = fDensityRatioNuD;
102 }
103 if(subType == fIonisation ||
104 subType == fPairProdByCharged ||
105 subType == fAnnihilation ||
106 subType == fAnnihilationToMuMu ||
107 subType == fAnnihilationToHadrons){
108 fProcessToDensity[processName] = fDensityRatioElD;
109 }
110 if(subType == fBremsstrahlung ||
111 subType == fNuclearStopping){
112 fProcessToDensity[processName] = fDensityRatioNuDElD;
113 }
114
115 if(subType == fCerenkov ||
116 subType == fScintillation ||
117 subType == fSynchrotronRadiation ||
118 subType == fTransitionRadiation){
119 fProcessToDensity[processName] = fDensityRatioNone;
120 }
121 if(subType == fRayleigh ||
122 subType == fPhotoElectricEffect ||
123 subType == fComptonScattering ||
124 subType == fGammaConversion ||
125 subType == fGammaConversionToMuMu){
126 fProcessToDensity[processName] = fDensityRatioNone;
127 }
128 break;
129 case fOptical:
130 fProcessToDensity[processName] = fDensityRatioNone;
131 break;
132 case fHadronic:
133 fProcessToDensity[processName] = fDensityRatioNuD;
134 break;
136 fProcessToDensity[processName] = fDensityRatioNuD;
137 break;
138 case fGeneral:
139 fProcessToDensity[processName] = fDensityRatioNone;
140 break;
141 case fDecay:
142 fProcessToDensity[processName] = fDensityRatioNone;
143 break;
145 fProcessToDensity[processName] = fDensityRatioNone;
146 break;
147 case fUserDefined:
148 fProcessToDensity[processName] = fDensityRatioNone;
149 break;
150 case fParallel:
151 fProcessToDensity[processName] = fDensityRatioNone;
152 break;
153 case fPhonon:
154 fProcessToDensity[processName] = fDensityRatioNone;
155 break;
156 case fUCN:
157 fProcessToDensity[processName] = fDensityRatioNone;
158 break;
159 default:
160 fProcessToDensity[processName] = fDensityRatioNone;
161 break;
162 }
163 }
164 }
165 fSetup = false;
166 }
167}
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
170
172G4ChannelingOptrChangeCrossSection::ProposeOccurenceBiasingOperation(const G4Track* track,
174 callingProcess)
175{
176 if ( track->GetDefinition() != fParticleToBias ) return 0;
177
178 G4double analogInteractionLength =
180 if ( analogInteractionLength > DBL_MAX/10. ) return 0;
181
182 G4double analogXS = 1./analogInteractionLength;
183
184 if(fChannelingID==-1){
185 fChannelingID = G4PhysicsModelCatalog::GetIndex("channeling");
186 }
187 G4ChannelingTrackData* trackdata =
188 (G4ChannelingTrackData*)(track->GetAuxiliaryTrackInformation(fChannelingID));
189 if(trackdata==nullptr) return 0;
190
191 G4double XStransformation = 1.;
192 auto search = fProcessToDensity.find(callingProcess->GetWrappedProcess()->GetProcessName());
193 if(search != fProcessToDensity.end()) {
194 switch (search->second) {
196 XStransformation = trackdata->GetDensity();
197 break;
198 case fDensityRatioNuD:
199 XStransformation = trackdata->GetNuD();
200 break;
201 case fDensityRatioElD:
202 XStransformation = trackdata->GetElD();
203 break;
205 return 0;
206 break;
208 return 0;
209 break;
210 default:
211 return 0;
212 break;
213 }
214 }
215 else{
216 XStransformation = trackdata->GetDensity();
217 }
218
219 G4BOptnChangeCrossSection* operation = fChangeCrossSectionOperations[callingProcess];
220 G4VBiasingOperation* previousOperation = callingProcess->GetPreviousOccurenceBiasingOperation();
221
222 if ( previousOperation == 0 ){
223 operation->SetBiasedCrossSection( XStransformation * analogXS );
224 operation->Sample();
225 }
226 else{
227 if ( previousOperation != operation ){
229 ed << " Logic problem in operation handling !" << G4endl;
230 G4Exception("G4ChannelingOptrChangeCrossSection::ProposeOccurenceBiasingOperation(...)",
231 "G4Channeling",
233 ed);
234 return 0;
235 }
236 if ( operation->GetInteractionOccured() ){
237 operation->SetBiasedCrossSection( XStransformation * analogXS );
238 operation->Sample();
239 }
240 else{
241 operation->UpdateForStep( callingProcess->GetPreviousStepSize() );
242 operation->SetBiasedCrossSection( XStransformation * analogXS );
243 operation->UpdateForStep( 0.0 );
244 }
245 }
246
247 return operation;
248
249}
250
251//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
252
253void G4ChannelingOptrChangeCrossSection::
254OperationApplied(const G4BiasingProcessInterface* callingProcess,
256 G4VBiasingOperation* occurenceOperationApplied,
257 G4double,
259 const G4VParticleChange* )
260{
261 G4BOptnChangeCrossSection* operation = fChangeCrossSectionOperations[callingProcess];
262 if ( operation == occurenceOperationApplied ) operation->SetInteractionOccured();
263}
264
265//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4BiasingAppliedCase
@ fGammaConversionToMuMu
@ fAnnihilationToHadrons
@ fBremsstrahlung
@ fCoulombScattering
@ fGammaConversion
@ fRayleigh
@ fIonisation
@ fPairProdByCharged
@ fSynchrotronRadiation
@ fCerenkov
@ fAnnihilationToMuMu
@ fScintillation
@ fNuclearStopping
@ fComptonScattering
@ fTransitionRadiation
@ fAnnihilation
@ fMultipleScattering
@ fPhotoElectricEffect
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ProcessType
@ fOptical
@ fPhonon
@ fParameterisation
@ fParallel
@ fUCN
@ fGeneral
@ fDecay
@ fElectromagnetic
@ fHadronic
@ fUserDefined
@ fTransportation
@ fPhotolepton_hadron
@ fNotDefined
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
void UpdateForStep(G4double stepLength)
void SetBiasedCrossSection(G4double xst, bool updateInteractionLength=false)
const G4BiasingProcessSharedData * GetSharedData() const
G4VBiasingOperation * GetPreviousOccurenceBiasingOperation() const
G4VProcess * GetWrappedProcess() const
const std::vector< const G4BiasingProcessInterface * > & GetPhysicsBiasingProcessInterfaces() const
G4ChannelingOptrChangeCrossSection(G4String particleToBias, G4String name="ChannelingChangeXS")
G4ProcessManager * GetProcessManager() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4int GetIndex(const G4String &)
G4ParticleDefinition * GetDefinition() const
G4VAuxiliaryTrackInformation * GetAuxiliaryTrackInformation(G4int idx) const
Definition: G4Track.cc:218
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:388
G4double GetCurrentInteractionLength() const
Definition: G4VProcess.hh:443
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
#define DBL_MAX
Definition: templates.hh:62