Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4VEmAdjointModel.hh
이 파일의 문서화 페이지로 가기
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: G4VEmAdjointModel.hh 100341 2016-10-18 08:02:25Z gcosmo $
27 //
29 // Module: G4VEMAdjointModel
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 // 10 September 2009 Move to a virtual class. L. Desorgher
40 // 1st April 2007 creation by L. Desorgher
41 //
42 //-------------------------------------------------------------
43 // Documentation:
44 // Base class for Adjoint EM model. It is based on the use of direct G4VEmModel.
45 //
46 
47 
48 #ifndef G4VEmAdjointModel_h
49 #define G4VEmAdjointModel_h 1
50 
51 #include "globals.hh"
52 #include "G4DynamicParticle.hh"
53 #include "G4ParticleDefinition.hh"
54 #include "G4MaterialCutsCouple.hh"
55 #include "G4Material.hh"
56 #include "G4Element.hh"
57 #include "G4ElementVector.hh"
58 #include "Randomize.hh"
59 #include "G4ParticleDefinition.hh"
60 #include "G4VEmModel.hh"
61 #include "G4Electron.hh"
62 #include "G4Gamma.hh"
63 #include "G4ProductionCutsTable.hh"
64 
65 class G4PhysicsTable;
66 class G4Region;
67 class G4VParticleChange;
68 class G4ParticleChange;
69 class G4Track;
70 class G4AdjointCSMatrix;
71 
73 {
74 
75 public: // public methods
76 
77  G4VEmAdjointModel(const G4String& nam);
78 
79  virtual ~G4VEmAdjointModel();
80 
81  //------------------------------------------------------------------------
82  // Virtual methods to be implemented for the sample secondaries concrete model
83  //------------------------------------------------------------------------
84 
85  //virtual void Initialise()=0;
86 
87  virtual void SampleSecondaries(const G4Track& aTrack,
88  G4bool IsScatProjToProjCase,
89  G4ParticleChange* fParticleChange)=0;
90 
91 
92  //------------------------------------------------------------------------
93  // Methods for adjoint processes; may be overwritten if needed;
94  //------------------------------------------------------------------------
95 
96 
97  virtual G4double AdjointCrossSection(const G4MaterialCutsCouple* aCouple,
98  G4double primEnergy,
99  G4bool IsScatProjToProjCase);
100 
102  G4double primEnergy,
103  G4bool IsScatProjToProjCase);
104 
106  G4double kinEnergyProj, // kinetic energy of the primary particle before the interaction
107  G4double kinEnergyProd, // kinetic energy of the secondary particle
108  G4double Z,
109  G4double A = 0.);
110 
112  G4double kinEnergyProj, // kinetic energy of the primary particle before the interaction
113  G4double kinEnergyScatProj, // kinetic energy of the primary particle after the interaction
114  G4double Z,
115  G4double A = 0.);
116 
117 
118 
120  const G4Material* aMaterial,
121  G4double kinEnergyProj, // kinetic energy of the primary particle before the interaction
122  G4double kinEnergyProd // kinetic energy of the secondary particle
123  );
124 
126  const G4Material* aMaterial,
127  G4double kinEnergyProj, // kinetic energy of the primary particle before the interaction
128  G4double kinEnergyScatProj // kinetic energy of the primary particle after the interaction
129  );
130 
131 
132  //Energy limits of adjoint secondary
133  //------------------
134 
139 
140 
141 
142  //Other Methods
143  //---------------
144 
145  void DefineCurrentMaterial(const G4MaterialCutsCouple* couple);
146 
147 
148  std::vector< std::vector< double>* > ComputeAdjointCrossSectionVectorPerAtomForSecond(
149  G4double kinEnergyProd,
150  G4double Z,
151  G4double A = 0.,
152  G4int nbin_pro_decade=10
153  );
154  std::vector< std::vector< double>* > ComputeAdjointCrossSectionVectorPerAtomForScatProj(
155  G4double kinEnergyProd,
156  G4double Z,
157  G4double A = 0.,
158  G4int nbin_pro_decade=10
159  );
160 
161  std::vector< std::vector< double>* > ComputeAdjointCrossSectionVectorPerVolumeForSecond(
162  G4Material* aMaterial,
163  G4double kinEnergyProd,
164  G4int nbin_pro_decade=10
165  );
166  std::vector< std::vector< double>* > ComputeAdjointCrossSectionVectorPerVolumeForScatProj(
167  G4Material* aMaterial,
168  G4double kinEnergyProd,
169  G4int nbin_pro_decade=10
170  );
171 
172 
173 
174  inline void SetCSMatrices(std::vector< G4AdjointCSMatrix* >* Vec1CSMatrix, std::vector< G4AdjointCSMatrix* >* Vec2CSMatrix){
177 
178 
179  };
180 
182 
184 
186 
188 
189  void SetHighEnergyLimit(G4double aVal);
190 
191  void SetLowEnergyLimit(G4double aVal);
192 
193  inline void DefineDirectEMModel(G4VEmModel* aModel){theDirectEMModel = aModel;}
194 
196 
199  }
200 
202 
204 
205  inline void SetUseMatrix(G4bool aBool) { UseMatrix = aBool;}
206 
207  inline void SetUseMatrixPerElement(G4bool aBool){ UseMatrixPerElement = aBool;}
209 
210  inline void SetApplyCutInRange(G4bool aBool){ ApplyCutInRange = aBool;}
211  inline G4bool GetUseMatrix() {return UseMatrix;}
215 
216  inline G4String GetName(){ return name;}
217  inline virtual void SetCSBiasingFactor(G4double aVal) {CS_biasing_factor = aVal;}
218 
221 
222 protected:
223 
224  //Some of them can be overriden by daughter classes
225 
226 
230 
231 
232 
233  //General methods to sample secondary energy
234  //--------------------------------------
235  G4double SampleAdjSecEnergyFromCSMatrix(size_t MatrixIndex,G4double prim_energy,G4bool IsScatProjToProjCase);
236  G4double SampleAdjSecEnergyFromCSMatrix(G4double prim_energy,G4bool IsScatProjToProjCase);
237  void SelectCSMatrix(G4bool IsScatProjToProjCase);
238 
239  virtual G4double SampleAdjSecEnergyFromDiffCrossSectionPerAtom(G4double prim_energy,G4bool IsScatProjToProjCase);
240 
241 
242 
243  //Post Step weight correction
244  //----------------------------
245  virtual void CorrectPostStepWeight(G4ParticleChange* fParticleChange,
246  G4double old_weight,
247  G4double adjointPrimKinEnergy,
248  G4double projectileKinEnergy,
249  G4bool IsScatProjToProjCase);
250 
251 
252 
253 
254 
255 
256 protected: //attributes
257 
260 
261 
262 
263 
264  //Name
265  //-----
266 
267  const G4String name;
268 
269  //Needed for CS integration at the initialisation phase
270  //-----------------------------------------------------
271 
278 
279  //for the adjoint simulation we need for each element or material:
280  //an adjoint CS Matrix
281  //-----------------------------
282 
283  std::vector< G4AdjointCSMatrix* >* pOnCSMatrixForProdToProjBackwardScattering;
284  std::vector< G4AdjointCSMatrix* >* pOnCSMatrixForScatProjToProjBackwardScattering;
285  std::vector<G4double> CS_Vs_ElementForScatProjToProjCase;
286  std::vector<G4double> CS_Vs_ElementForProdToProjCase;
287 
291 
292  //particle definition
293  //------------------
294 
299 
300  //Prestep energy
301  //-------------
303 
304  //Current couple material
305  //----------------------
313 
314  //For ions
315  //---------
318 
319  //Energy limits
320  //-------------
321 
324 
325  //Cross Section biasing factor
326  //---------------------------
328 
329  //Type of Model with Matrix or not
330  //--------------------------------
332  G4bool UseMatrixPerElement; //other possibility is per Material
334 
335  //Index of Cross section matrices to be used
336  //------------
338 
339  size_t model_index;
340 
341  //This is needed for the forced interaction where part of the weight correction
342  // is given outside the model while the secondary are created in the model
343  //The weight should be fixed before adding the secondary
346 
347 };
348 
349 
350 #endif
351 
G4double DiffCrossSectionFunction1(G4double kinEnergyProj)
virtual G4double GetSecondAdjEnergyMinForProdToProjCase(G4double PrimAdjEnergy)
G4VEmAdjointModel(const G4String &nam)
G4double lastAdjointCSForScatProjToProjCase
G4bool GetSecondPartOfSameType()
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerAtomForSecond(G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
G4ParticleDefinition * GetAdjointEquivalentOfDirectSecondaryParticleDefinition()
G4double GetLowEnergyLimit()
G4double DiffCrossSectionFunction2(G4double kinEnergyProj)
virtual G4double GetSecondAdjEnergyMinForScatProjToProjCase(G4double PrimAdjEnergy, G4double Tcut=0)
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerAtomForScatProj(G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
void DefineDirectEMModel(G4VEmModel *aModel)
void SetCSMatrices(std::vector< G4AdjointCSMatrix * > *Vec1CSMatrix, std::vector< G4AdjointCSMatrix * > *Vec2CSMatrix)
void SetApplyCutInRange(G4bool aBool)
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerVolumeForScatProj(G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
G4Material * currentMaterial
void SetUseOnlyOneMatrixForAllElements(G4bool aBool)
G4double GetHighEnergyLimit()
G4bool UseOnlyOneMatrixForAllElements
G4double lastAdjointCSForProdToProjCase
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerVolumeForSecond(G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
G4double additional_weight_correction_factor_for_post_step_outside_model
G4double currentTcutForDirectSecond
G4MaterialCutsCouple * currentCouple
G4double kinEnergyProjForIntegration
Float_t Z
G4double kinEnergyProdForIntegration
G4ParticleDefinition * theAdjEquivOfDirectSecondPartDef
std::vector< G4AdjointCSMatrix * > * pOnCSMatrixForProdToProjBackwardScattering
void SetAdjointEquivalentOfDirectSecondaryParticleDefinition(G4ParticleDefinition *aPart)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void SetCorrectWeightForPostStepInModel(G4bool aBool)
G4ParticleDefinition * theAdjEquivOfDirectPrimPartDef
virtual G4double DiffCrossSectionPerAtomPrimToSecond(G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
virtual G4double DiffCrossSectionPerVolumePrimToSecond(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
void SelectCSMatrix(G4bool IsScatProjToProjCase)
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase(G4double PrimAdjEnergy)
G4bool GetUseOnlyOneMatrixForAllElements()
double A(double temperature)
size_t indexOfUsedCrossSectionMatrix
G4bool correct_weight_for_post_step_in_model
virtual void CorrectPostStepWeight(G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool IsScatProjToProjCase)
G4VEmModel * theDirectEMModel
virtual G4double DiffCrossSectionPerVolumePrimToScatPrim(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyScatProj)
void SetSecondPartOfSameType(G4bool aBool)
int G4int
Definition: G4Types.hh:78
std::vector< G4double > CS_Vs_ElementForProdToProjCase
G4double DiffCrossSectionPerVolumeFunctionForIntegrationOverEkinProj(G4double EkinProd)
G4double currentTcutForDirectPrim
void SetAdditionalWeightCorrectionFactorForPostStepOutsideModel(G4double factor)
virtual G4double SampleAdjSecEnergyFromDiffCrossSectionPerAtom(G4double prim_energy, G4bool IsScatProjToProjCase)
virtual void SetCSBiasingFactor(G4double aVal)
void SetLowEnergyLimit(G4double aVal)
G4double kinEnergyScatProjForIntegration
virtual G4double DiffCrossSectionPerAtomPrimToScatPrim(G4double kinEnergyProj, G4double kinEnergyScatProj, G4double Z, G4double A=0.)
G4double SampleAdjSecEnergyFromCSMatrix(size_t MatrixIndex, G4double prim_energy, G4bool IsScatProjToProjCase)
void SetHighEnergyLimit(G4double aVal)
G4ParticleDefinition * theDirectPrimaryPartDef
G4Material * SelectedMaterial
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
void SetUseMatrixPerElement(G4bool aBool)
virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase(G4double PrimAdjEnergy)
G4VParticleChange * pParticleChange
G4ParticleDefinition * GetAdjointEquivalentOfDirectPrimaryParticleDefinition()
virtual G4double GetAdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
virtual void SampleSecondaries(const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)=0
std::vector< G4double > CS_Vs_ElementForScatProjToProjCase
std::vector< G4AdjointCSMatrix * > * pOnCSMatrixForScatProjToProjBackwardScattering
void SetAdjointEquivalentOfDirectPrimaryParticleDefinition(G4ParticleDefinition *aPart)
void SetUseMatrix(G4bool aBool)
void DefineCurrentMaterial(const G4MaterialCutsCouple *couple)
virtual ~G4VEmAdjointModel()