Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointAlongStepWeightCorrection.hh
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//
27/////////////////////////////////////////////////////////////////////////////////
28// Class: G4AdjointAlongStepWeightCorrection
29// Author: L. Desorgher
30// Organisation: SpaceIT GmbH
31// Contract: ESA contract 21435/08/NL/AT
32// Customer: ESA/ESTEC
33/////////////////////////////////////////////////////////////////////////////////
34//
35// CHANGE HISTORY
36// --------------
37// ChangeHistory:
38// 10 May 2007 creation by L. Desorgher
39// October 2009 implementation of the mode where the total adjoint and forward cross sections are equivalent. L. Desorgher
40//
41//-------------------------------------------------------------
42// Documentation:
43// Continuous processes acting on adjoint particles to correct continuously their weight during the adjoint reverse tracking.
44// Thi process is needed whene the adjoint cross section are not scaled such that the total adjoint cross section match the total forward cross section.
45// By default the mode where the total adjoint cross section is equal to the total forward cross section is used an therefore this along step weight
46// correction factor is 1.
47// However in some cases (some energy ranges) the total forward cross section or the total adjoint cross section can be null, in this case the along step
48// weight correction is neede and is given by exp(-(Sigma_tot_adj-Sigma_tot_fwd).dx)
49//
50//
51//
52
53
54#ifndef G4AdjointAlongStepWeightCorrection_h
55#define G4AdjointAlongStepWeightCorrection_h 1
56
58#include "globals.hh"
59#include "G4Material.hh"
61#include "G4Track.hh"
62#include "G4ParticleChange.hh"
63
64class G4Step;
66
67
68
70{
71public:
72
73 G4AdjointAlongStepWeightCorrection(const G4String& name = "ContinuousWeightCorrection",
75
77
78
79protected:
80 virtual G4double GetContinuousStepLimit(const G4Track& track,
81 G4double previousStepSize,
82 G4double currentMinimumStep,
83 G4double& currentSafety);
84
85
86 //------------------------------------------------------------------------
87 // Generic methods common to all processes
88 //------------------------------------------------------------------------
89public:
90
92
94
96
97
98private:
99
100 void DefineMaterial(const G4MaterialCutsCouple* couple);
101
102
103
106
107
108
109protected:
110
112
113
114private:
115
116 const G4Material* currentMaterial;
117 const G4MaterialCutsCouple* currentCouple;
118 size_t currentMaterialIndex;
119 G4double preStepKinEnergy;
120
121};
122
123inline void G4AdjointAlongStepWeightCorrection::DefineMaterial(
124 const G4MaterialCutsCouple* couple)
125{
126 if(couple != currentCouple) {
127 currentCouple = couple;
128 currentMaterial = couple->GetMaterial();
129 currentMaterialIndex = couple->GetIndex();
130
131 }
132}
133
134
135
136#endif
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
void BuildPhysicsTable(const G4ParticleDefinition &)
void PreparePhysicsTable(const G4ParticleDefinition &)
virtual G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
const G4Material * GetMaterial() const
Definition: G4Step.hh:62