source: trunk/source/processes/electromagnetic/adjoint/include/G4VEmAdjointModel.hh @ 966

Last change on this file since 966 was 966, checked in by garnier, 15 years ago

fichier ajoutes

File size: 11.8 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/////////////////////////////////////////////////////////////////////////////////
27//      Module:         G4VEMAdjointModel.hh
28//      Author:         L. Desorgher
29//      Date:           1st April 2007
30//      Organisation:   SpaceIT GmbH
31//      Customer:       ESA/ESTEC
32/////////////////////////////////////////////////////////////////////////////////
33//
34// CHANGE HISTORY
35// --------------
36//      ChangeHistory:
37//              1st April 2007 creation by L. Desorgher                 
38//
39//-------------------------------------------------------------
40//      Documentation:
41//              Base class for Adjoint model
42//
43
44#ifndef G4VEmAdjointModel_h
45#define G4VEmAdjointModel_h 1
46
47#include "globals.hh"
48#include "G4DynamicParticle.hh"
49#include "G4ParticleDefinition.hh"
50#include "G4MaterialCutsCouple.hh"
51#include "G4Material.hh"
52#include "G4Element.hh"
53#include "G4ElementVector.hh"
54#include "Randomize.hh"
55#include "G4ParticleDefinition.hh"
56#include "G4VEmModel.hh"
57#include "G4Electron.hh"
58#include "G4Gamma.hh"
59#include "G4ProductionCutsTable.hh"
60
61class G4PhysicsTable;
62class G4Region;
63class G4VParticleChange;
64class G4ParticleChange;
65class G4Track;
66class G4AdjointCSMatrix;
67
68class G4VEmAdjointModel
69{
70
71public:
72
73  G4VEmAdjointModel(const G4String& nam);
74
75  virtual ~G4VEmAdjointModel();
76
77  //------------------------------------------------------------------------
78  // Virtual methods to be implemented for the concrete model
79  //------------------------------------------------------------------------
80
81  //virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&) = 0;
82
83
84  virtual void SampleSecondaries(const G4Track& aTrack,
85                                G4bool IsScatProjToProjCase,
86                                G4ParticleChange* fParticleChange);
87 
88 
89
90  //------------------------------------------------------------------------
91  // Methods for adjoint processes; may be overwritten if needed; 
92  //------------------------------------------------------------------------
93  virtual G4double AdjointCrossSection(const G4MaterialCutsCouple* aCouple,
94                                             G4double primEnergy,
95                                             G4bool IsScatProjToProjCase);
96                               
97  virtual G4double DiffCrossSectionPerAtomPrimToSecond(
98                                      G4double kinEnergyProj,  // kinetic energy of the primary particle before the interaction
99                                      G4double kinEnergyProd, // kinetic energy of the secondary particle
100                                      G4double Z, 
101                                      G4double A = 0.);
102                                     
103  virtual G4double DiffCrossSectionPerAtomPrimToScatPrim( 
104                                      G4double kinEnergyProj,  // kinetic energy of the primary particle before the interaction
105                                      G4double kinEnergyScatProj, // kinetic energy of the primary particle after the interaction
106                                      G4double Z, 
107                                      G4double A = 0.);
108 
109 
110 
111  virtual G4double DiffCrossSectionPerVolumePrimToSecond(
112                                      const G4Material* aMaterial,
113                                      G4double kinEnergyProj,  // kinetic energy of the primary particle before the interaction
114                                      G4double kinEnergyProd // kinetic energy of the secondary particle
115                                      );
116                                     
117  virtual G4double DiffCrossSectionPerVolumePrimToScatPrim(
118                                      const G4Material* aMaterial, 
119                                      G4double kinEnergyProj,  // kinetic energy of the primary particle before the interaction
120                                      G4double kinEnergyScatProj // kinetic energy of the primary particle after the interaction
121                                      );
122 
123                                     
124 
125 
126  G4double DiffCrossSectionFunction1(G4double kinEnergyProj);
127  G4double DiffCrossSectionMoller(G4double kinEnergyProj,G4double kinEnergyProd);
128 
129  G4double DiffCrossSectionFunction2(G4double kinEnergyProj);
130 
131  std::vector< std::vector< G4double >* >  ComputeAdjointCrossSectionVectorPerAtomForSecond(     
132                                G4double kinEnergyProd,
133                                G4double Z, 
134                                G4double A = 0.,
135                                G4int nbin_pro_decade=10
136                                );
137  std::vector< std::vector< G4double >* >  ComputeAdjointCrossSectionVectorPerAtomForScatProj(     
138                                G4double kinEnergyProd,
139                                G4double Z, 
140                                G4double A = 0.,
141                                G4int nbin_pro_decade=10
142                                );
143 
144  std::vector< std::vector< G4double >* >  ComputeAdjointCrossSectionVectorPerVolumeForSecond(     
145                                G4Material* aMaterial,
146                                G4double kinEnergyProd,
147                                G4int nbin_pro_decade=10
148                                );
149  std::vector< std::vector< G4double >* >  ComputeAdjointCrossSectionVectorPerVolumeForScatProj(     
150                                G4Material* aMaterial,
151                                G4double kinEnergyProd,
152                                G4int nbin_pro_decade=10
153                                );
154                               
155 
156                                                       
157  virtual G4double SampleAdjSecEnergyFromCSMatrix(size_t MatrixIndex,G4double prim_energy,G4bool IsScatProjToProjCase);     
158  virtual G4double SampleAdjSecEnergyFromDiffCrossSectionPerAtom(G4double prim_energy,G4bool IsScatProjToProjCase);     
159  void CorrectPostStepWeight(G4ParticleChange* fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy);     
160 
161  //Set/Get methods
162  //------------------
163 
164  virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase(G4double PrimAdjEnergy);
165  virtual G4double GetSecondAdjEnergyMinForScatProjToProjCase(G4double PrimAdjEnergy,G4double Tcut=0);
166  virtual G4double GetSecondAdjEnergyMaxForProdToProjCase(G4double PrimAdjEnergy);
167  virtual G4double GetSecondAdjEnergyMinForProdToProjCase(G4double PrimAdjEnergy);
168 
169  virtual void SetCSBiasingFactor(G4double aVal) {CS_biasing_factor = aVal;}
170 
171 
172 
173 
174 
175 
176public: 
177 
178  inline void SetCSMatrices(std::vector< G4AdjointCSMatrix* >* Vec1CSMatrix, std::vector< G4AdjointCSMatrix* >* Vec2CSMatrix){
179                                 pOnCSMatrixForProdToProjBackwardScattering = Vec1CSMatrix;
180                                 pOnCSMatrixForScatProjToProjBackwardScattering = Vec2CSMatrix;
181                                 
182       
183  };
184 
185  inline G4ParticleDefinition* GetAdjointEquivalentOfDirectPrimaryParticleDefinition(){return theAdjEquivOfDirectPrimPartDef;}
186 
187  inline G4ParticleDefinition* GetAdjointEquivalentOfDirectSecondaryParticleDefinition(){return theAdjEquivOfDirectSecondPartDef;}     
188 
189  inline G4double GetHighEnergyLimit(){return HighEnergyLimit;}
190 
191  inline G4double GetLowEnergyLimit(){return LowEnergyLimit;}
192 
193  inline void SetHighEnergyLimit(G4double aVal){HighEnergyLimit=aVal;}
194 
195  inline void SetLowEnergyLimit(G4double aVal){LowEnergyLimit=aVal;}
196 
197  inline void SetCorrectWeightMode(G4bool aBool){CorrectWeightMode=aBool;};
198 
199  inline void SetApplyBiasing(G4bool aBool){ApplyBiasing=aBool;};
200 
201                             
202  inline void DefineDirectEMModel(G4VEmModel* aModel){theDirectEMModel = aModel;}
203 
204  inline void SetAdjointEquivalentOfDirectPrimaryParticleDefinition(G4ParticleDefinition* aPart){
205        theAdjEquivOfDirectPrimPartDef=aPart;
206        if (theAdjEquivOfDirectPrimPartDef->GetParticleName() =="adj_e-")
207                                        theDirectPrimaryPartDef=G4Electron::Electron();
208        if (theAdjEquivOfDirectPrimPartDef->GetParticleName() =="adj_gamma")
209                                        theDirectPrimaryPartDef=G4Gamma::Gamma();
210 
211  }
212 
213  inline void SetAdjointEquivalentOfDirectSecondaryParticleDefinition(G4ParticleDefinition* aPart){
214        theAdjEquivOfDirectSecondPartDef =aPart;
215  }
216 
217  inline void SetSecondPartOfSameType(G4bool aBool){second_part_of_same_type =aBool;}
218 
219  bool GetSecondPartOfSameType(){return second_part_of_same_type;}
220 
221  inline void SetUseMatrix(G4bool aBool) { UseMatrix = aBool;}
222 
223  inline void SetUseMatrixPerElement(G4bool aBool){ UseMatrixPerElement = aBool;}
224  inline void SetUseOnlyOneMatrixForAllElements(G4bool aBool){ UseOnlyOneMatrixForAllElements = aBool;}
225 
226  inline void SetApplyCutInRange(G4bool aBool){ ApplyCutInRange = aBool;} 
227  inline void SetIsIonisation(G4bool aBool){ IsIonisation = aBool;} 
228 
229  inline G4bool GetUseMatrix() {return UseMatrix;}
230  inline G4bool GetUseMatrixPerElement(){ return UseMatrixPerElement;} 
231  inline G4bool GetUseOnlyOneMatrixForAllElements(){ return UseOnlyOneMatrixForAllElements;} 
232  inline G4bool GetApplyCutInRange(){ return ApplyCutInRange;} 
233  void  DefineCurrentMaterial(const G4MaterialCutsCouple* couple);
234 
235  inline G4String GetName(){ return name;} 
236
237private: //Methods
238 
239 
240 
241                                   
242protected:
243 
244  G4VEmModel* theDirectEMModel;
245  G4VParticleChange*  pParticleChange;
246 
247
248protected:
249
250  //  hide assignment operator
251  G4VEmAdjointModel & operator=(const  G4VEmAdjointModel &right);
252  G4VEmAdjointModel(const  G4VEmAdjointModel&);
253 
254 
255  //Name
256  //-----
257 
258  const G4String  name;
259 
260  //Needed for CS integration at the initialisation phase
261  //-----------------------------------------------------
262 
263  G4int ASelectedNucleus;
264  G4int ZSelectedNucleus;
265  G4Material* SelectedMaterial;
266  G4double kinEnergyProdForIntegration;
267  G4double kinEnergyScatProjForIntegration;
268 
269 
270  //for the adjoint simulation  we need for each element or material:
271  //an adjoint CS Matrix
272  //-----------------------------
273 
274  std::vector< G4AdjointCSMatrix* >* pOnCSMatrixForProdToProjBackwardScattering;
275  std::vector< G4AdjointCSMatrix* >* pOnCSMatrixForScatProjToProjBackwardScattering;
276  std::vector<double> CS_Vs_ElementForScatProjToProjCase;
277  std::vector<double> CS_Vs_ElementForProdToProjCase;
278 
279  G4double lastCS;
280 
281  //particle definition
282  //------------------
283 
284  G4ParticleDefinition* theAdjEquivOfDirectPrimPartDef;
285  G4ParticleDefinition* theAdjEquivOfDirectSecondPartDef;
286  G4ParticleDefinition* theDirectPrimaryPartDef;
287  G4bool second_part_of_same_type;
288 
289  //Current couple material
290  //----------------------
291  G4Material*  currentMaterial;
292  G4MaterialCutsCouple* currentCouple;
293  size_t   currentMaterialIndex; 
294  size_t   currentCoupleIndex; 
295  G4double currentTcutForDirectPrim;
296  G4double currentTcutForDirectSecond;
297  G4bool ApplyCutInRange;
298 
299 
300  //CorrectWeightMode
301  //------------------
302 
303  bool CorrectWeightMode;
304 
305  //Apply biasing
306  //------------
307 
308  bool ApplyBiasing;
309 
310 
311 
312 
313 
314  //Energy limits
315  //-------------
316 
317  G4double HighEnergyLimit;
318  G4double LowEnergyLimit; 
319 
320 
321  //Cross Section biasing factor
322  //---------------------------
323  G4double CS_biasing_factor;
324 
325 
326  //Type of Model with Matrix or not
327  //--------------------------------
328   bool UseMatrix;
329   bool UseMatrixPerElement; //other possibility is per Material
330   bool UseOnlyOneMatrixForAllElements;
331   bool IsIonisation;
332   
333   
334   
335
336 
337};
338
339
340#endif
341
Note: See TracBrowser for help on using the repository browser.