Geant4-11
G4VEmAdjointModel.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//
27// Class: G4VEMAdjointModel
28// Author: L. Desorgher
29// Organisation: SpaceIT GmbH
30//
31// Base class for Adjoint EM model. It is based on the use of direct
32// G4VEmModel.
34
35#ifndef G4VEmAdjointModel_h
36#define G4VEmAdjointModel_h 1
37
38#include "globals.hh"
40#include "G4VEmModel.hh"
41
44class G4Material;
47class G4Region;
48class G4Track;
49
51{
52 public:
53 explicit G4VEmAdjointModel(const G4String& nam);
54
55 virtual ~G4VEmAdjointModel();
56
57 //------------------------------------------------------------------------
58 // Virtual methods to be implemented for the sample secondaries concrete model
59 //------------------------------------------------------------------------
60
61 virtual void SampleSecondaries(const G4Track& aTrack, G4bool isScatProjToProj,
62 G4ParticleChange* fParticleChange) = 0;
63
64 //------------------------------------------------------------------------
65 // Methods for adjoint processes
66 //------------------------------------------------------------------------
67
69 G4double primEnergy,
70 G4bool isScatProjToProj);
71
72 // The implementation of the DiffCrossSection... here are correct for
73 // energy loss process. For the photoelectric and Compton scattering
74 // the method should be redefined
76 G4double kinEnergyProj, // kin energy of primary before interaction
77 G4double kinEnergyProd, // kinetic energy of the secondary particle
78 G4double Z, G4double A = 0.);
79
81 G4double kinEnergyProj, // kin energy of primary before interaction
82 G4double kinEnergyScatProj, // kin energy of primary after interaction
83 G4double Z, G4double A = 0.);
84
86 const G4Material* aMaterial,
87 G4double kinEnergyProj, // kin energy of primary before interaction
88 G4double kinEnergyProd // kinetic energy of secondary particle
89 );
90
92 const G4Material* aMaterial,
93 G4double kinEnergyProj, // kin energy of primary before interaction
94 G4double kinEnergyScatProj // kinetic energy of primary after interaction
95 );
96
97 // Energy limits of adjoint secondary
98 //------------------
99
101 G4double primAdjEnergy);
102
104 G4double primAdjEnergy, G4double tcut = 0.);
105
107
109
110 // Other Methods
111 //---------------
112
114
115 std::vector<std::vector<double>*>
117 G4double Z, G4double A = 0.,
118 G4int nbin_pro_decade = 10);
119
120 std::vector<std::vector<double>*>
122 G4double kinEnergyProd, G4double Z, G4double A = 0.,
123 G4int nbin_pro_decade = 10);
124
125 std::vector<std::vector<double>*>
127 G4Material* aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade = 10);
128
129 std::vector<std::vector<double>*>
131 G4Material* aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade = 10);
132
133 inline void SetCSMatrices(std::vector<G4AdjointCSMatrix*>* Vec1CSMatrix,
134 std::vector<G4AdjointCSMatrix*>* Vec2CSMatrix)
135 {
136 fCSMatrixProdToProjBackScat = Vec1CSMatrix;
137 fCSMatrixProjToProjBackScat = Vec2CSMatrix;
138 };
139
142 {
144 }
145
148 {
150 }
151
153
155
156 void SetHighEnergyLimit(G4double aVal);
157
158 void SetLowEnergyLimit(G4double aVal);
159
160 inline void DefineDirectEMModel(G4VEmModel* aModel) { fDirectModel = aModel; }
161
163 G4ParticleDefinition* aPart);
164
167 {
169 }
170
172 {
173 fSecondPartSameType = aBool;
174 }
175
177
178 inline void SetUseMatrix(G4bool aBool) { fUseMatrix = aBool; }
179
181 {
182 fUseMatrixPerElement = aBool;
183 }
184
186 {
188 }
189
190 inline void SetApplyCutInRange(G4bool aBool) { fApplyCutInRange = aBool; }
191
192 inline G4bool GetUseMatrix() { return fUseMatrix; }
193
195
197 {
199 }
200
202
203 inline G4String GetName() { return fName; }
204
205 inline virtual void SetCSBiasingFactor(G4double aVal)
206 {
207 fCsBiasingFactor = aVal;
208 }
209
211 {
212 fInModelWeightCorr = aBool;
213 }
214
216 G4double factor)
217 {
218 fOutsideWeightFactor = factor;
219 }
220
223
224 protected:
226
228
229 // General methods to sample secondary energy
230 G4double SampleAdjSecEnergyFromCSMatrix(size_t MatrixIndex,
231 G4double prim_energy,
232 G4bool isScatProjToProj);
233
235 G4bool isScatProjToProj);
236
237 void SelectCSMatrix(G4bool isScatProjToProj);
238
240 G4double prim_energy, G4bool isScatProjToProj);
241
242 // Post Step weight correction
243 virtual void CorrectPostStepWeight(G4ParticleChange* fParticleChange,
244 G4double old_weight,
245 G4double adjointPrimKinEnergy,
246 G4double projectileKinEnergy,
247 G4bool isScatProjToProj);
248
251
253
257
258 // particle definition
262
263 // adjoint CS matrix for each element or material
264 std::vector<G4AdjointCSMatrix*>* fCSMatrixProdToProjBackScat = nullptr;
265 std::vector<G4AdjointCSMatrix*>* fCSMatrixProjToProjBackScat = nullptr;
266
267 std::vector<G4double> fElementCSScatProjToProj;
268 std::vector<G4double> fElementCSProdToProj;
269
272
276
278
281
282 // Energy limits
285
286 // Cross Section biasing factor
288
289 // [1] This is needed for the forced interaction where part of the weight
290 // correction is given outside the model while the secondary are created in
291 // the model. The weight should be fixed before adding the secondary
293
294 // Needed for CS integration at the initialisation phase
297
298 size_t fCSMatrixUsed = 0; // Index of crosssection matrices used
299
302 false; // correct_weight_for_post_step_in_model, see [1]
303
305
306 // Type of Model with Matrix or not
308 G4bool fUseMatrixPerElement = false; // other possibility is per Material
310};
311
312#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
virtual G4double GetSecondAdjEnergyMaxForProdToProj(G4double primAdjEnergy)
G4VEmAdjointModel(G4VEmAdjointModel &)=delete
void SetUseMatrixPerElement(G4bool aBool)
G4double fLastAdjointCSForScatProjToProj
std::vector< G4AdjointCSMatrix * > * fCSMatrixProdToProjBackScat
virtual void SampleSecondaries(const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange)=0
G4ParticleDefinition * fAdjEquivDirectSecondPart
virtual G4double DiffCrossSectionPerVolumePrimToScatPrim(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyScatProj)
void SetAdjointEquivalentOfDirectSecondaryParticleDefinition(G4ParticleDefinition *aPart)
G4bool GetSecondPartOfSameType()
void SetCSMatrices(std::vector< G4AdjointCSMatrix * > *Vec1CSMatrix, std::vector< G4AdjointCSMatrix * > *Vec2CSMatrix)
virtual ~G4VEmAdjointModel()
G4double fKinEnergyScatProjForIntegration
G4double DiffCrossSectionFunction1(G4double kinEnergyProj)
void SetUseMatrix(G4bool aBool)
G4bool GetUseOnlyOneMatrixForAllElements()
virtual G4double GetSecondAdjEnergyMinForScatProjToProj(G4double primAdjEnergy, G4double tcut=0.)
G4double fKinEnergyProdForIntegration
G4ParticleDefinition * fAdjEquivDirectPrimPart
G4MaterialCutsCouple * fCurrentCouple
virtual G4double GetSecondAdjEnergyMaxForScatProjToProj(G4double primAdjEnergy)
G4Material * fCurrentMaterial
std::vector< G4AdjointCSMatrix * > * fCSMatrixProjToProjBackScat
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerAtomForSecond(G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
virtual void SetCSBiasingFactor(G4double aVal)
virtual G4double DiffCrossSectionPerVolumePrimToSecond(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
virtual void CorrectPostStepWeight(G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool isScatProjToProj)
std::vector< G4double > fElementCSScatProjToProj
G4ParticleDefinition * GetAdjointEquivalentOfDirectPrimaryParticleDefinition()
G4double GetLowEnergyLimit()
G4double SampleAdjSecEnergyFromCSMatrix(size_t MatrixIndex, G4double prim_energy, G4bool isScatProjToProj)
virtual G4double DiffCrossSectionPerAtomPrimToScatPrim(G4double kinEnergyProj, G4double kinEnergyScatProj, G4double Z, G4double A=0.)
G4AdjointCSManager * fCSManager
void SetLowEnergyLimit(G4double aVal)
void SetAdditionalWeightCorrectionFactorForPostStepOutsideModel(G4double factor)
void DefineCurrentMaterial(const G4MaterialCutsCouple *couple)
G4ParticleDefinition * GetAdjointEquivalentOfDirectSecondaryParticleDefinition()
void SetSecondPartOfSameType(G4bool aBool)
void SelectCSMatrix(G4bool isScatProjToProj)
const G4String fName
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool isScatProjToProj)
G4VEmModel * fDirectModel
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerVolumeForSecond(G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
G4Material * fSelectedMaterial
void SetCorrectWeightForPostStepInModel(G4bool aBool)
void SetUseOnlyOneMatrixForAllElements(G4bool aBool)
virtual G4double GetSecondAdjEnergyMinForProdToProj(G4double primAdjEnergy)
G4double DiffCrossSectionFunction2(G4double kinEnergyProj)
G4ParticleDefinition * fDirectPrimaryPart
void SetAdjointEquivalentOfDirectPrimaryParticleDefinition(G4ParticleDefinition *aPart)
G4VEmAdjointModel(const G4String &nam)
void DefineDirectEMModel(G4VEmModel *aModel)
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerVolumeForScatProj(G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
virtual G4double DiffCrossSectionPerAtomPrimToSecond(G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
std::vector< G4double > fElementCSProdToProj
void SetHighEnergyLimit(G4double aVal)
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerAtomForScatProj(G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
G4double GetHighEnergyLimit()
G4double fLastAdjointCSForProdToProj
virtual G4double SampleAdjSecEnergyFromDiffCrossSectionPerAtom(G4double prim_energy, G4bool isScatProjToProj)
void SetApplyCutInRange(G4bool aBool)
G4VEmAdjointModel & operator=(const G4VEmAdjointModel &right)=delete