Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointCSManager.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: G4AdjointCSManager
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// 1st April 2007 creation by L. Desorgher
39//
40// September-October 2009. Implementation of the mode where the adjoint cross sections are scaled such that the total used adjoint cross sections is in
41// most of the cases equal to the total forward cross section. L.Desorgher
42//
43//-------------------------------------------------------------
44// Documentation:
45// Is responsible for the management of all adjoint cross sections matrices, and for the computation of the total forward and adjoint cross sections.
46// Total adjoint and forward cross sections are needed to correct the weight of a particle after a tracking step or after the occurrence of a reverse reaction.
47// It is also used to sample an adjoint secondary from a given adjoint cross section matrix.
48//
49#ifndef G4AdjointCSManager_h
50#define G4AdjointCSManager_h 1
51
52#include"globals.hh"
53#include<vector>
56
59class G4Material;
61class G4Element;
62class G4VEmProcess;
64class G4PhysicsTable;
65
66////////////////////////////////////////////////////////////////////////////////
67//
69{
70
72
73 public:
76
77 public:
79
80 //Registration of the different models and processes
81
83
84 void RegisterEmProcess(G4VEmProcess* aProcess, G4ParticleDefinition* aPartDef);
85
87
89
90 //Building of the CS Matrices and Total Forward and Adjoint LambdaTables
91 //----------------------------------------------------------------------
92
95
96
97 //Get TotalCrossSections form Total Lambda Tables, Needed for Weight correction and scaling of the
98 //-------------------------------------------------
100 const G4MaterialCutsCouple* aCouple);
102 const G4MaterialCutsCouple* aCouple);
103
104 G4double GetAdjointSigma(G4double Ekin_nuc, size_t index_model,G4bool is_scat_proj_to_proj,
105 const G4MaterialCutsCouple* aCouple);
106
108 const G4MaterialCutsCouple* aCouple, G4double& emin_adj, G4double& emin_fwd);
110 const G4MaterialCutsCouple* aCouple, G4double& e_sigma_max, G4double& sigma_max);
112 const G4MaterialCutsCouple* aCouple, G4double& e_sigma_max, G4double& sigma_max);
113
114
115
116 //CrossSection Correction 1 or FwdCS/AdjCS following the G4boolean value of forward_CS_is_used and forward_CS_mode
117 //-------------------------------------------------
118 G4double GetCrossSectionCorrection(G4ParticleDefinition* aPartDef,G4double PreStepEkin,const G4MaterialCutsCouple* aCouple, G4bool& fwd_is_used, G4double& fwd_TotCS);
119
120
121 //Cross section mode
122 //------------------
123 inline void SetFwdCrossSectionMode(G4bool aBool){forward_CS_mode=aBool;}
124
125
126 //Weight correction
127 //------------------
128
130 const G4MaterialCutsCouple* aCouple, G4double step_length);
132
133
134
135
136 //Method Called by the adjoint model to get there CS, if not precised otherwise
137 //-------------------------------
138
140 G4VEmAdjointModel* aModel,
141 G4double PrimEnergy,
142 G4double Tcut,
143 G4bool IsScatProjToProjCase,
144 std::vector<G4double>&
145 AdjointCS_for_each_element);
146
147 //Method Called by the adjoint model to sample the secondary energy form the CS matrix
148 //--------------------------------------------------------------------------------
150 G4VEmAdjointModel* aModel,
151 G4double PrimEnergy,
152 G4double Tcut,
153 G4bool IsScatProjToProjCase);
154
155
156 //Total Adjoint CS is computed at initialisation phase
157 //-----------------------------------------------------
159
160
161
162
165
166 //inline
167 inline void SetTmin(G4double aVal){Tmin=aVal;}
168 inline void SetTmax(G4double aVal){Tmax=aVal;}
169 inline void SetNbins(G4int aInt){nbins=aInt;}
170 inline void SetIon(G4ParticleDefinition* adjIon,
171 G4ParticleDefinition* fwdIon) {theAdjIon=adjIon; theFwdIon =fwdIon;}
172
173
174 private:
175 static G4ThreadLocal G4AdjointCSManager* theInstance;
176 std::vector< std::vector<G4AdjointCSMatrix*> > theAdjointCSMatricesForScatProjToProj; //x dim is for G4VAdjointEM*, y dim is for elements
177 std::vector< std::vector<G4AdjointCSMatrix*> > theAdjointCSMatricesForProdToProj;
178 std::vector< G4VEmAdjointModel*> listOfAdjointEMModel;
179
180 std::vector<G4AdjointCSMatrix*>
181 BuildCrossSectionsMatricesForAGivenModelAndElement(G4VEmAdjointModel* aModel,
182 G4int Z,
183 G4int A,
184 G4int nbin_pro_decade);
185
186 std::vector<G4AdjointCSMatrix*>
187 BuildCrossSectionsMatricesForAGivenModelAndMaterial(G4VEmAdjointModel* aModel,
188 G4Material* aMaterial,
189 G4int nbin_pro_decade);
190
191
192 G4Material* lastMaterial;
193 G4double lastPrimaryEnergy;
194 G4double lastTcut;
195 std::vector< size_t> listOfIndexOfAdjointEMModelInAction;
196 std::vector< G4bool> listOfIsScatProjToProjCase;
197 std::vector< std::vector<G4double> > lastAdjointCSVsModelsAndElements;
198 G4bool CrossSectionMatrixesAreBuilt;
199 size_t currentParticleIndex;
200 G4ParticleDefinition* currentParticleDef;
201
202 //total adjoint and total forward cross section table in function of material and in function of adjoint particle type
203 //--------------------------------------------------------------------------------------------------------------------
204 std::vector<G4PhysicsTable*> theTotalForwardSigmaTableVector;
205 std::vector<G4PhysicsTable*> theTotalAdjointSigmaTableVector;
206 std::vector< std::vector<G4double> > EminForFwdSigmaTables;
207 std::vector< std::vector<G4double> > EminForAdjSigmaTables;
208 std::vector< std::vector<G4double> > EkinofFwdSigmaMax;
209 std::vector< std::vector<G4double> > EkinofAdjSigmaMax;
210 G4bool TotalSigmaTableAreBuilt;
211
212 //Sigma tavle for each G4VAdjointEMModel
213 std::vector<G4PhysicsTable*> listSigmaTableForAdjointModelScatProjToProj;
214 std::vector<G4PhysicsTable*> listSigmaTableForAdjointModelProdToProj;
215
216 //list of forward G4VEMLossProcess and of G4VEMProcess for the different adjoint particle
217 //--------------------------------------------------------------
218 std::vector< std::vector<G4VEmProcess*>* > listOfForwardEmProcess;
219 std::vector< std::vector<G4VEnergyLossProcess*>* > listOfForwardEnergyLossProcess;
220
221 //list of adjoint particles considered
222 //--------------------------------------------------------------
223 std::vector< G4ParticleDefinition*> theListOfAdjointParticlesInAction;
224
225 G4double Tmin,Tmax;
226 G4int nbins;
227
228 //Current material
229 //----------------
230 G4MaterialCutsCouple* currentCouple;
231 G4Material* currentMaterial;
232 size_t currentMatIndex;
233
234 G4int verbose;
235
236 //Two CS mode are possible :forward_CS_mode = false the Adjoint CS are used as it is implying a AlongStep Weight Correction.
237 // :forward_CS_mode = true the Adjoint CS are scaled to have the total adjoint CS eual to the fwd one implying a PostStep Weight Correction.
238 // For energy range where the total FwdCS or the total adjoint CS are null, the scaling is not possble and
239 // forward_CS_is_used is set to false
240 //--------------------------------------------
241 G4bool forward_CS_is_used;
242 G4bool forward_CS_mode;
243
244 //Adj and Fwd CS values for re-use
245 //------------------------
246
247 G4double PreadjCS,PostadjCS;
248 G4double PrefwdCS,PostfwdCS;
249 G4double LastEkinForCS;
250 G4double LastCSCorrectionFactor;
251 G4ParticleDefinition* lastPartDefForCS;
252
253 //Ion
254 //----------------
255 G4ParticleDefinition* theAdjIon; //at the moment Only one ion can be considered by simulation
256 G4ParticleDefinition* theFwdIon;
257 G4double massRatio;
258
259 private:
261 void DefineCurrentMaterial(const G4MaterialCutsCouple* couple);
262 void DefineCurrentParticle(const G4ParticleDefinition* aPartDef);
263 G4double ComputeAdjointCS(G4double aPrimEnergy, G4AdjointCSMatrix* anAdjointCSMatrix, G4double Tcut);
264 size_t eindex;
265
266};
267#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
size_t RegisterEmAdjointModel(G4VEmAdjointModel *)
void RegisterAdjointParticle(G4ParticleDefinition *aPartDef)
G4ParticleDefinition * GetForwardParticleEquivalent(G4ParticleDefinition *theAdjPartDef)
void GetMaxAdjTotalCS(G4ParticleDefinition *aPartDef, const G4MaterialCutsCouple *aCouple, G4double &e_sigma_max, G4double &sigma_max)
G4double GetCrossSectionCorrection(G4ParticleDefinition *aPartDef, G4double PreStepEkin, const G4MaterialCutsCouple *aCouple, G4bool &fwd_is_used, G4double &fwd_TotCS)
G4double GetTotalForwardCS(G4ParticleDefinition *aPartDef, G4double Ekin, const G4MaterialCutsCouple *aCouple)
void SetIon(G4ParticleDefinition *adjIon, G4ParticleDefinition *fwdIon)
G4double ComputeTotalAdjointCS(const G4MaterialCutsCouple *aMatCutCouple, G4ParticleDefinition *aPart, G4double PrimEnergy)
G4ParticleDefinition * GetAdjointParticleEquivalent(G4ParticleDefinition *theFwdPartDef)
G4Element * SampleElementFromCSMatrices(G4Material *aMaterial, G4VEmAdjointModel *aModel, G4double PrimEnergy, G4double Tcut, G4bool IsScatProjToProjCase)
void SetNbins(G4int aInt)
G4double GetContinuousWeightCorrection(G4ParticleDefinition *aPartDef, G4double PreStepEkin, G4double AfterStepEkin, const G4MaterialCutsCouple *aCouple, G4double step_length)
void SetTmin(G4double aVal)
G4double GetPostStepWeightCorrection()
void GetEminForTotalCS(G4ParticleDefinition *aPartDef, const G4MaterialCutsCouple *aCouple, G4double &emin_adj, G4double &emin_fwd)
static G4AdjointCSManager * GetAdjointCSManager()
G4double GetTotalAdjointCS(G4ParticleDefinition *aPartDef, G4double Ekin, const G4MaterialCutsCouple *aCouple)
G4double GetAdjointSigma(G4double Ekin_nuc, size_t index_model, G4bool is_scat_proj_to_proj, const G4MaterialCutsCouple *aCouple)
void RegisterEmProcess(G4VEmProcess *aProcess, G4ParticleDefinition *aPartDef)
void SetTmax(G4double aVal)
void RegisterEnergyLossProcess(G4VEnergyLossProcess *aProcess, G4ParticleDefinition *aPartDef)
void GetMaxFwdTotalCS(G4ParticleDefinition *aPartDef, const G4MaterialCutsCouple *aCouple, G4double &e_sigma_max, G4double &sigma_max)
G4double ComputeAdjointCS(G4Material *aMaterial, G4VEmAdjointModel *aModel, G4double PrimEnergy, G4double Tcut, G4bool IsScatProjToProjCase, std::vector< G4double > &AdjointCS_for_each_element)
void SetFwdCrossSectionMode(G4bool aBool)
#define G4ThreadLocal
Definition: tls.hh:77