source: trunk/source/processes/electromagnetic/adjoint/include/G4AdjointCSManager.hh@ 1197

Last change on this file since 1197 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 10.7 KB
Line 
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// $Id: G4AdjointCSManager.hh,v 1.4 2009/11/20 10:31:20 ldesorgh Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29/////////////////////////////////////////////////////////////////////////////////
30// Class: G4AdjointCSManager
31// Author: L. Desorgher
32// Organisation: SpaceIT GmbH
33// Contract: ESA contract 21435/08/NL/AT
34// Customer: ESA/ESTEC
35/////////////////////////////////////////////////////////////////////////////////
36//
37// CHANGE HISTORY
38// --------------
39// ChangeHistory:
40// 1st April 2007 creation by L. Desorgher
41//
42// September-October 2009. Implementation of the mode where the adjoint cross sections are scaled such that the total used adjoint cross sections is in
43// most of the cases equal to the total forward cross section. L.Desorgher
44//
45//-------------------------------------------------------------
46// Documentation:
47// Is responsible for the management of all adjoint cross sections matrices, and for the computation of the total forward and adjoint cross sections.
48// Total adjoint and forward cross sections are needed to correct the weight of a particle after a tracking step or after the occurence of a reverse reaction.
49// It is also used to sample an adjoint secondary from a given adjoint cross section matrix.
50//
51#ifndef G4AdjointCSManager_h
52#define G4AdjointCSManager_h 1
53
54#include"globals.hh"
55#include<vector>
56#include"G4AdjointCSMatrix.hh"
57
58
59class G4VEmAdjointModel;
60class G4MaterialCutsCouple;
61class G4Material;
62class G4ParticleDefinition;
63class G4Element;
64class G4VEmProcess;
65class G4VEnergyLossProcess;
66class G4PhysicsTable;
67
68////////////////////////////////////////////////////////////////////////////////
69//
70class G4AdjointCSManager
71{
72
73 public:
74 ~G4AdjointCSManager();
75 static G4AdjointCSManager* GetAdjointCSManager();
76
77 public:
78 G4int GetNbProcesses();
79
80 //Registration of the different models and processes
81
82 void RegisterEmAdjointModel(G4VEmAdjointModel*);
83
84 void RegisterEmProcess(G4VEmProcess* aProcess, G4ParticleDefinition* aPartDef);
85
86 void RegisterEnergyLossProcess(G4VEnergyLossProcess* aProcess, G4ParticleDefinition* aPartDef);
87
88 void RegisterAdjointParticle(G4ParticleDefinition* aPartDef);
89
90 //Building of the CS Matrices and Total Forward and Adjoint LambdaTables
91 //----------------------------------------------------------------------
92
93 void BuildCrossSectionMatrices();
94 void BuildTotalSigmaTables();
95
96
97 //Get TotalCrossSections form Total Lambda Tables, Needed for Weight correction and scaling of the
98 //-------------------------------------------------
99 G4double GetTotalAdjointCS(G4ParticleDefinition* aPartDef, G4double Ekin,
100 const G4MaterialCutsCouple* aCouple);
101 G4double GetTotalForwardCS(G4ParticleDefinition* aPartDef, G4double Ekin,
102 const G4MaterialCutsCouple* aCouple);
103
104
105 void GetEminForTotalCS(G4ParticleDefinition* aPartDef,
106 const G4MaterialCutsCouple* aCouple, G4double& emin_adj, G4double& emin_fwd);
107 void GetMaxFwdTotalCS(G4ParticleDefinition* aPartDef,
108 const G4MaterialCutsCouple* aCouple, G4double& e_sigma_max, G4double& sigma_max);
109 void GetMaxAdjTotalCS(G4ParticleDefinition* aPartDef,
110 const G4MaterialCutsCouple* aCouple, G4double& e_sigma_max, G4double& sigma_max);
111
112
113
114 //CrossSection Correction 1 or FwdCS/AdjCS following the G4boolean value of forward_CS_is_used and forward_CS_mode
115 //-------------------------------------------------
116 G4double GetCrossSectionCorrection(G4ParticleDefinition* aPartDef,G4double PreStepEkin,const G4MaterialCutsCouple* aCouple, G4bool& fwd_is_used, G4double& fwd_TotCS);
117
118
119 //Cross section mode
120 //------------------
121 inline void SetFwdCrossSectionMode(G4bool aBool){forward_CS_mode=aBool;}
122
123
124 //Weight correction
125 //------------------
126
127 G4double GetContinuousWeightCorrection(G4ParticleDefinition* aPartDef, G4double PreStepEkin,G4double AfterStepEkin,
128 const G4MaterialCutsCouple* aCouple, G4double step_length);
129 G4double GetPostStepWeightCorrection();
130
131
132
133
134 //Method Called by the adjoint model to get there CS, if not precised otherwise
135 //-------------------------------
136
137 G4double ComputeAdjointCS(G4Material* aMaterial,
138 G4VEmAdjointModel* aModel,
139 G4double PrimEnergy,
140 G4double Tcut,
141 G4bool IsScatProjToProjCase,
142 std::vector<G4double>&
143 AdjointCS_for_each_element);
144
145 //Method Called by the adjoint model to sample the secondary energy form the CS matrix
146 //--------------------------------------------------------------------------------
147 G4Element* SampleElementFromCSMatrices(G4Material* aMaterial,
148 G4VEmAdjointModel* aModel,
149 G4double PrimEnergy,
150 G4double Tcut,
151 G4bool IsScatProjToProjCase);
152
153
154 //Total Adjoint CS is computed at initialisation phase
155 //-----------------------------------------------------
156 G4double ComputeTotalAdjointCS(const G4MaterialCutsCouple* aMatCutCouple,G4ParticleDefinition* aPart,G4double PrimEnergy);
157
158
159
160
161 G4ParticleDefinition* GetAdjointParticleEquivalent(G4ParticleDefinition* theFwdPartDef);
162 G4ParticleDefinition* GetForwardParticleEquivalent(G4ParticleDefinition* theAdjPartDef);
163
164 //inline
165 inline void SetTmin(G4double aVal){Tmin=aVal;}
166 inline void SetTmax(G4double aVal){Tmax=aVal;}
167 inline void SetNbins(G4int aInt){nbins=aInt;}
168 inline void SetIon(G4ParticleDefinition* adjIon,
169 G4ParticleDefinition* fwdIon) {theAdjIon=adjIon; theFwdIon =fwdIon;}
170
171
172 private:
173 static G4AdjointCSManager* theInstance;
174 std::vector< std::vector<G4AdjointCSMatrix*> > theAdjointCSMatricesForScatProjToProj; //x dim is for G4VAdjointEM* while y dim is for elements
175 std::vector< std::vector<G4AdjointCSMatrix*> > theAdjointCSMatricesForProdToProj;
176 std::vector< G4VEmAdjointModel*> listOfAdjointEMModel;
177
178 std::vector<G4AdjointCSMatrix*>
179 BuildCrossSectionsMatricesForAGivenModelAndElement(G4VEmAdjointModel* aModel,
180 G4int Z,
181 G4int A,
182 G4int nbin_pro_decade);
183
184 std::vector<G4AdjointCSMatrix*>
185 BuildCrossSectionsMatricesForAGivenModelAndMaterial(G4VEmAdjointModel* aModel,
186 G4Material* aMaterial,
187 G4int nbin_pro_decade);
188
189
190 G4Material* lastMaterial;
191 G4double lastPrimaryEnergy;
192 G4double lastTcut;
193 std::vector< size_t> listOfIndexOfAdjointEMModelInAction;
194 std::vector< G4bool> listOfIsScatProjToProjCase;
195 std::vector< std::vector<G4double> > lastAdjointCSVsModelsAndElements;
196 G4bool CrossSectionMatrixesAreBuilt;
197 size_t currentParticleIndex;
198 G4ParticleDefinition* currentParticleDef;
199
200 //total adjoint and total forward cross section table in function of material and in function of adjoint particle type
201 //--------------------------------------------------------------------------------------------------------------------
202 std::vector<G4PhysicsTable*> theTotalForwardSigmaTableVector;
203 std::vector<G4PhysicsTable*> theTotalAdjointSigmaTableVector;
204 std::vector< std::vector<G4double> > EminForFwdSigmaTables;
205 std::vector< std::vector<G4double> > EminForAdjSigmaTables;
206 std::vector< std::vector<G4double> > EkinofFwdSigmaMax;
207 std::vector< std::vector<G4double> > EkinofAdjSigmaMax;
208
209
210
211
212 //list of forward G4VEMLossProcess and of G4VEMProcess for the different adjoint particle
213 //--------------------------------------------------------------
214 std::vector< std::vector<G4VEmProcess*>* > listOfForwardEmProcess;
215 std::vector< std::vector<G4VEnergyLossProcess*>* > listOfForwardEnergyLossProcess;
216
217 //list of adjoint particles considered
218 //--------------------------------------------------------------
219 std::vector< G4ParticleDefinition*> theListOfAdjointParticlesInAction;
220
221
222 G4double Tmin,Tmax;
223 G4int nbins;
224
225
226 //Current material
227 //----------------
228 G4MaterialCutsCouple* currentCouple;
229 G4Material* currentMaterial;
230 size_t currentMatIndex;
231
232 G4int verbose;
233
234
235
236
237 //Two CS mode are possible :forward_CS_mode = false the Adjoint CS are used as it is implying a AlongStep Weight Correction.
238 // :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.
239 // For energy range where the total FwdCS or the total adjoint CS are null, the scaling is not possble and
240 // forward_CS_is_used is set to false
241 //--------------------------------------------
242 G4bool forward_CS_is_used;
243 G4bool forward_CS_mode;
244
245 //Adj and Fwd CS values for re-use
246 //------------------------
247
248 G4double PreadjCS,PostadjCS;
249 G4double PrefwdCS,PostfwdCS;
250 G4double LastEkinForCS;
251 G4double LastCSCorrectionFactor;
252 G4ParticleDefinition* lastPartDefForCS;
253
254 //Ion
255 //----------------
256 G4ParticleDefinition* theAdjIon; //at the moment Only one ion can be considered by simulation
257 G4ParticleDefinition* theFwdIon;
258 G4double massRatio;
259
260
261
262
263 private:
264 G4AdjointCSManager();
265 void DefineCurrentMaterial(const G4MaterialCutsCouple* couple);
266 void DefineCurrentParticle(const G4ParticleDefinition* aPartDef);
267 G4double ComputeAdjointCS(G4double aPrimEnergy, G4AdjointCSMatrix* anAdjointCSMatrix, G4double Tcut);
268
269};
270#endif
Note: See TracBrowser for help on using the repository browser.