Geant4  10.02.p02
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 // $Id: G4AdjointCSManager.hh 93569 2015-10-26 14:53:21Z gcosmo $
27 //
29 // Class: G4AdjointCSManager
30 // Author: L. Desorgher
31 // Organisation: SpaceIT GmbH
32 // Contract: ESA contract 21435/08/NL/AT
33 // Customer: ESA/ESTEC
35 //
36 // CHANGE HISTORY
37 // --------------
38 // ChangeHistory:
39 // 1st April 2007 creation by L. Desorgher
40 //
41 // September-October 2009. Implementation of the mode where the adjoint cross sections are scaled such that the total used adjoint cross sections is in
42 // most of the cases equal to the total forward cross section. L.Desorgher
43 //
44 //-------------------------------------------------------------
45 // Documentation:
46 // Is responsible for the management of all adjoint cross sections matrices, and for the computation of the total forward and adjoint cross sections.
47 // 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.
48 // It is also used to sample an adjoint secondary from a given adjoint cross section matrix.
49 //
50 #ifndef G4AdjointCSManager_h
51 #define G4AdjointCSManager_h 1
52 
53 #include"globals.hh"
54 #include<vector>
55 #include"G4AdjointCSMatrix.hh"
57 
58 class G4VEmAdjointModel;
60 class G4Material;
62 class G4Element;
63 class G4VEmProcess;
65 class G4PhysicsTable;
66 
68 //
70 {
71 
73 
74  public:
77 
78  public:
80 
81  //Registration of the different models and processes
82 
84 
85  void RegisterEmProcess(G4VEmProcess* aProcess, G4ParticleDefinition* aPartDef);
86 
88 
90 
91  //Building of the CS Matrices and Total Forward and Adjoint LambdaTables
92  //----------------------------------------------------------------------
93 
95  void BuildTotalSigmaTables();
96 
97 
98  //Get TotalCrossSections form Total Lambda Tables, Needed for Weight correction and scaling of the
99  //-------------------------------------------------
101  const G4MaterialCutsCouple* aCouple);
103  const G4MaterialCutsCouple* aCouple);
104 
105  G4double GetAdjointSigma(G4double Ekin_nuc, size_t index_model,G4bool is_scat_proj_to_proj,
106  const G4MaterialCutsCouple* aCouple);
107 
109  const G4MaterialCutsCouple* aCouple, G4double& emin_adj, G4double& emin_fwd);
110  void GetMaxFwdTotalCS(G4ParticleDefinition* aPartDef,
111  const G4MaterialCutsCouple* aCouple, G4double& e_sigma_max, G4double& sigma_max);
112  void GetMaxAdjTotalCS(G4ParticleDefinition* aPartDef,
113  const G4MaterialCutsCouple* aCouple, G4double& e_sigma_max, G4double& sigma_max);
114 
115 
116 
117  //CrossSection Correction 1 or FwdCS/AdjCS following the G4boolean value of forward_CS_is_used and forward_CS_mode
118  //-------------------------------------------------
119  G4double GetCrossSectionCorrection(G4ParticleDefinition* aPartDef,G4double PreStepEkin,const G4MaterialCutsCouple* aCouple, G4bool& fwd_is_used, G4double& fwd_TotCS);
120 
121 
122  //Cross section mode
123  //------------------
124  inline void SetFwdCrossSectionMode(G4bool aBool){forward_CS_mode=aBool;}
125 
126 
127  //Weight correction
128  //------------------
129 
130  G4double GetContinuousWeightCorrection(G4ParticleDefinition* aPartDef, G4double PreStepEkin,G4double AfterStepEkin,
131  const G4MaterialCutsCouple* aCouple, G4double step_length);
133 
134 
135 
136 
137  //Method Called by the adjoint model to get there CS, if not precised otherwise
138  //-------------------------------
139 
141  G4VEmAdjointModel* aModel,
142  G4double PrimEnergy,
143  G4double Tcut,
144  G4bool IsScatProjToProjCase,
145  std::vector<G4double>&
146  AdjointCS_for_each_element);
147 
148  //Method Called by the adjoint model to sample the secondary energy form the CS matrix
149  //--------------------------------------------------------------------------------
151  G4VEmAdjointModel* aModel,
152  G4double PrimEnergy,
153  G4double Tcut,
154  G4bool IsScatProjToProjCase);
155 
156 
157  //Total Adjoint CS is computed at initialisation phase
158  //-----------------------------------------------------
159  G4double ComputeTotalAdjointCS(const G4MaterialCutsCouple* aMatCutCouple,G4ParticleDefinition* aPart,G4double PrimEnergy);
160 
161 
162 
163 
166 
167  //inline
168  inline void SetTmin(G4double aVal){Tmin=aVal;}
169  inline void SetTmax(G4double aVal){Tmax=aVal;}
170  inline void SetNbins(G4int aInt){nbins=aInt;}
171  inline void SetIon(G4ParticleDefinition* adjIon,
172  G4ParticleDefinition* fwdIon) {theAdjIon=adjIon; theFwdIon =fwdIon;}
173 
174 
175  private:
177  std::vector< std::vector<G4AdjointCSMatrix*> > theAdjointCSMatricesForScatProjToProj; //x dim is for G4VAdjointEM*, y dim is for elements
178  std::vector< std::vector<G4AdjointCSMatrix*> > theAdjointCSMatricesForProdToProj;
179  std::vector< G4VEmAdjointModel*> listOfAdjointEMModel;
180 
181  std::vector<G4AdjointCSMatrix*>
183  G4int Z,
184  G4int A,
185  G4int nbin_pro_decade);
186 
187  std::vector<G4AdjointCSMatrix*>
189  G4Material* aMaterial,
190  G4int nbin_pro_decade);
191 
192 
197  std::vector< G4bool> listOfIsScatProjToProjCase;
198  std::vector< std::vector<G4double> > lastAdjointCSVsModelsAndElements;
202 
203  //total adjoint and total forward cross section table in function of material and in function of adjoint particle type
204  //--------------------------------------------------------------------------------------------------------------------
205  std::vector<G4PhysicsTable*> theTotalForwardSigmaTableVector;
206  std::vector<G4PhysicsTable*> theTotalAdjointSigmaTableVector;
207  std::vector< std::vector<G4double> > EminForFwdSigmaTables;
208  std::vector< std::vector<G4double> > EminForAdjSigmaTables;
209  std::vector< std::vector<G4double> > EkinofFwdSigmaMax;
210  std::vector< std::vector<G4double> > EkinofAdjSigmaMax;
212 
213  //Sigma tavle for each G4VAdjointEMModel
214  std::vector<G4PhysicsTable*> listSigmaTableForAdjointModelScatProjToProj;
215  std::vector<G4PhysicsTable*> listSigmaTableForAdjointModelProdToProj;
216 
217  //list of forward G4VEMLossProcess and of G4VEMProcess for the different adjoint particle
218  //--------------------------------------------------------------
219  std::vector< std::vector<G4VEmProcess*>* > listOfForwardEmProcess;
220  std::vector< std::vector<G4VEnergyLossProcess*>* > listOfForwardEnergyLossProcess;
221 
222  //list of adjoint particles considered
223  //--------------------------------------------------------------
224  std::vector< G4ParticleDefinition*> theListOfAdjointParticlesInAction;
225 
228 
229  //Current material
230  //----------------
234 
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  //--------------------------------------------
244 
245  //Adj and Fwd CS values for re-use
246  //------------------------
247 
253 
254  //Ion
255  //----------------
256  G4ParticleDefinition* theAdjIon; //at the moment Only one ion can be considered by simulation
259 
260  private:
262  void DefineCurrentMaterial(const G4MaterialCutsCouple* couple);
263  void DefineCurrentParticle(const G4ParticleDefinition* aPartDef);
264  G4double ComputeAdjointCS(G4double aPrimEnergy, G4AdjointCSMatrix* anAdjointCSMatrix, G4double Tcut);
265  size_t eindex;
266 
267 };
268 #endif
G4double GetTotalAdjointCS(G4ParticleDefinition *aPartDef, G4double Ekin, const G4MaterialCutsCouple *aCouple)
G4double GetCrossSectionCorrection(G4ParticleDefinition *aPartDef, G4double PreStepEkin, const G4MaterialCutsCouple *aCouple, G4bool &fwd_is_used, G4double &fwd_TotCS)
void SetTmin(G4double aVal)
G4double GetAdjointSigma(G4double Ekin_nuc, size_t index_model, G4bool is_scat_proj_to_proj, const G4MaterialCutsCouple *aCouple)
G4ParticleDefinition * theFwdIon
G4double ComputeTotalAdjointCS(const G4MaterialCutsCouple *aMatCutCouple, G4ParticleDefinition *aPart, G4double PrimEnergy)
G4double GetPostStepWeightCorrection()
std::vector< G4PhysicsTable * > theTotalAdjointSigmaTableVector
G4MaterialCutsCouple * currentCouple
std::vector< G4PhysicsTable * > listSigmaTableForAdjointModelProdToProj
void GetMaxFwdTotalCS(G4ParticleDefinition *aPartDef, const G4MaterialCutsCouple *aCouple, G4double &e_sigma_max, G4double &sigma_max)
size_t RegisterEmAdjointModel(G4VEmAdjointModel *)
void DefineCurrentMaterial(const G4MaterialCutsCouple *couple)
void SetIon(G4ParticleDefinition *adjIon, G4ParticleDefinition *fwdIon)
void RegisterEnergyLossProcess(G4VEnergyLossProcess *aProcess, G4ParticleDefinition *aPartDef)
#define G4ThreadLocal
Definition: tls.hh:89
std::vector< size_t > listOfIndexOfAdjointEMModelInAction
void SetFwdCrossSectionMode(G4bool aBool)
int G4int
Definition: G4Types.hh:78
std::vector< G4PhysicsTable * > listSigmaTableForAdjointModelScatProjToProj
std::vector< G4PhysicsTable * > theTotalForwardSigmaTableVector
void RegisterAdjointParticle(G4ParticleDefinition *aPartDef)
std::vector< std::vector< G4AdjointCSMatrix * > > theAdjointCSMatricesForProdToProj
std::vector< std::vector< G4double > > EminForFwdSigmaTables
static G4ThreadLocal G4AdjointCSManager * theInstance
double A(double temperature)
G4ParticleDefinition * GetAdjointParticleEquivalent(G4ParticleDefinition *theFwdPartDef)
std::vector< std::vector< G4double > > EkinofAdjSigmaMax
std::vector< std::vector< G4AdjointCSMatrix * > > theAdjointCSMatricesForScatProjToProj
std::vector< G4ParticleDefinition * > theListOfAdjointParticlesInAction
bool G4bool
Definition: G4Types.hh:79
G4ParticleDefinition * theAdjIon
std::vector< std::vector< G4VEnergyLossProcess * > * > listOfForwardEnergyLossProcess
void SetNbins(G4int aInt)
void SetTmax(G4double aVal)
std::vector< G4AdjointCSMatrix * > BuildCrossSectionsMatricesForAGivenModelAndMaterial(G4VEmAdjointModel *aModel, G4Material *aMaterial, G4int nbin_pro_decade)
G4double GetTotalForwardCS(G4ParticleDefinition *aPartDef, G4double Ekin, const G4MaterialCutsCouple *aCouple)
std::vector< std::vector< G4double > > EkinofFwdSigmaMax
std::vector< std::vector< G4double > > EminForAdjSigmaTables
void GetEminForTotalCS(G4ParticleDefinition *aPartDef, const G4MaterialCutsCouple *aCouple, G4double &emin_adj, G4double &emin_fwd)
G4ParticleDefinition * GetForwardParticleEquivalent(G4ParticleDefinition *theAdjPartDef)
G4double GetContinuousWeightCorrection(G4ParticleDefinition *aPartDef, G4double PreStepEkin, G4double AfterStepEkin, const G4MaterialCutsCouple *aCouple, G4double step_length)
G4ParticleDefinition * lastPartDefForCS
void RegisterEmProcess(G4VEmProcess *aProcess, G4ParticleDefinition *aPartDef)
std::vector< G4bool > listOfIsScatProjToProjCase
std::vector< std::vector< G4VEmProcess * > * > listOfForwardEmProcess
void GetMaxAdjTotalCS(G4ParticleDefinition *aPartDef, const G4MaterialCutsCouple *aCouple, G4double &e_sigma_max, G4double &sigma_max)
G4Material * currentMaterial
G4ParticleDefinition * currentParticleDef
std::vector< G4VEmAdjointModel * > listOfAdjointEMModel
std::vector< std::vector< G4double > > lastAdjointCSVsModelsAndElements
double G4double
Definition: G4Types.hh:76
G4Element * SampleElementFromCSMatrices(G4Material *aMaterial, G4VEmAdjointModel *aModel, G4double PrimEnergy, G4double Tcut, G4bool IsScatProjToProjCase)
void DefineCurrentParticle(const G4ParticleDefinition *aPartDef)
static G4AdjointCSManager * GetAdjointCSManager()
std::vector< G4AdjointCSMatrix * > BuildCrossSectionsMatricesForAGivenModelAndElement(G4VEmAdjointModel *aModel, G4int Z, G4int A, G4int nbin_pro_decade)
G4double ComputeAdjointCS(G4Material *aMaterial, G4VEmAdjointModel *aModel, G4double PrimEnergy, G4double Tcut, G4bool IsScatProjToProjCase, std::vector< G4double > &AdjointCS_for_each_element)