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

Last change on this file since 1358 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

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