source: trunk/source/processes/electromagnetic/lowenergy/include/G4PenelopeIonisationModel.hh @ 1064

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

maj sur la beta de geant 4.9.3

File size: 6.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: G4PenelopeIonisationModel.hh,v 1.1 2008/12/04 14:12:09 pandola Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-00 $
28//
29// Author: Luciano Pandola
30//
31// History:
32// -----------
33// 26 Nov 2008   L. Pandola   1st implementation. Migration from EM process
34//                            to EM model. Physics is unchanged.
35//
36// -------------------------------------------------------------------
37//
38// Class description:
39// Low Energy Electromagnetic Physics, e+ and e- ionisation
40// with Penelope Model
41// -------------------------------------------------------------------
42
43#ifndef G4PENELOPEIONISATIONMODEL_HH
44#define G4PENELOPEIONISATIONMODEL_HH 1
45
46#include "globals.hh"
47#include "G4VEmModel.hh"
48#include "G4DataVector.hh"
49#include "G4ParticleChangeForLoss.hh"
50#include "G4VCrossSectionHandler.hh"
51#include "G4PhysicsLogVector.hh"
52#include "G4AtomicDeexcitation.hh"
53
54class G4ParticleDefinition;
55class G4DynamicParticle;
56class G4MaterialCutsCouple;
57class G4Material;
58class G4VEMDataSet;
59
60class G4PenelopeIonisationModel : public G4VEmModel
61{
62
63public:
64 
65  G4PenelopeIonisationModel(const G4ParticleDefinition* p=0,
66                         const G4String& processName ="PenelopeIoni");
67 
68  virtual ~G4PenelopeIonisationModel();
69
70  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
71
72  virtual G4double CrossSectionPerVolume(const G4Material* material,
73                                         const G4ParticleDefinition* theParticle,
74                                         G4double kineticEnergy,
75                                         G4double cutEnergy,
76                                         G4double maxEnergy = DBL_MAX);
77                                         
78  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
79                                 const G4MaterialCutsCouple*,
80                                 const G4DynamicParticle*,
81                                 G4double tmin,
82                                 G4double maxEnergy);
83                                   
84  virtual G4double ComputeDEDXPerVolume(const G4Material*,
85                               const G4ParticleDefinition*,
86                               G4double kineticEnergy,
87                               G4double cutEnergy);
88                                 
89
90  void SetUseAtomicDeexcitation(G4bool value){fUseAtomicDeexcitation = value;};         
91  G4bool GetUseAtomicDeexcitation(){return fUseAtomicDeexcitation;};
92
93  void SetVerbosityLevel(G4int lev){verboseLevel = lev;};
94  G4int GetVerbosityLevel(){return verboseLevel;};
95
96
97protected:
98  G4ParticleChangeForLoss* fParticleChange;
99
100private:
101 
102  G4PenelopeIonisationModel & operator=(const G4PenelopeIonisationModel &right);
103  G4PenelopeIonisationModel(const G4PenelopeIonisationModel&);
104
105 
106  //Intrinsic energy limits of the model: cannot be extended by the parent process
107  G4double fIntrinsicLowEnergyLimit;
108  G4double fIntrinsicHighEnergyLimit;
109
110  G4bool fUseAtomicDeexcitation;
111
112  G4int verboseLevel;
113
114  G4bool isInitialised;
115 
116  G4double CalculateDeltaFermi(G4double kinEnergy ,G4int Z,
117                               G4double electronVolumeDensity);
118       
119  //Methods and variables to calculate final state
120  void CalculateDiscreteForElectrons(G4double kinEnergy,G4double cutoffEnergy,
121                                     G4int Z,G4double electronVolumeDensity);
122  void CalculateDiscreteForPositrons(G4double kinEnergy,G4double cutoffEnergy,
123                                     G4int Z,G4double electronVolumeDensity);
124
125  G4AtomicDeexcitation deexcitationManager;
126  G4double kineticEnergy1;
127  G4double cosThetaPrimary;
128  G4double energySecondary;
129  G4double cosThetaSecondary;
130  G4int iOsc;                             
131
132  //These methods are used to calculate the hard-cross section (namely they
133  //return the hard/total cross section)
134  G4double CalculateCrossSectionsRatio(G4double kinEnergy,
135                                       G4double cutoffEnergy,
136                                       G4int Z, 
137                                       G4double electronVolumeDensity,
138                                       const G4ParticleDefinition*);
139  //In fact the total cross section (hard+soft) is read from file
140  //The following methods give the cross section contribution (hard and soft) from each
141  //individual oscillator
142  std::pair<G4double,G4double> CrossSectionsRatioForElectrons(G4double kineticEnergy,
143                                                              G4double resEnergy,
144                                                              G4double densityCorrection,
145                                                              G4double cutoffEnergy);
146
147  std::pair<G4double,G4double> CrossSectionsRatioForPositrons(G4double kineticEnergy,
148                                                              G4double resEnergy,
149                                                              G4double densityCorrection,
150                                                              G4double cutoffEnergy);
151 
152  G4VCrossSectionHandler* crossSectionHandler;
153 
154  //These methods are used to calculate the stopping power up to the cutoff
155  //for each individual oscillator
156  G4double ComputeStoppingPowerForElectrons(G4double kinEnergy,
157                                            G4double cutEnergy,
158                                            G4double deltaFermi,
159                                            G4double resEnergy);
160
161  G4double ComputeStoppingPowerForPositrons(G4double kinEnergy,
162                                            G4double cutEnergy,
163                                            G4double deltaFermi,
164                                            G4double resEnergy);
165 
166 
167  //Parameters of atomic shells
168  void ReadData();
169  std::map<G4int,G4DataVector*> *ionizationEnergy;
170  std::map<G4int,G4DataVector*> *resonanceEnergy;
171  std::map<G4int,G4DataVector*> *occupationNumber;
172  std::map<G4int,G4DataVector*> *shellFlag;
173 
174  //Mean free path table. This will become obsolete! For now I need something to store
175  //cross sections and to sample a random atom
176  std::vector<G4VEMDataSet*>* theXSTable;
177  std::vector<G4VEMDataSet*>* BuildCrossSectionTable(const G4ParticleDefinition*);
178  G4int SampleRandomAtom(const G4MaterialCutsCouple*,G4double energy) const;
179
180};
181
182#endif
183
Note: See TracBrowser for help on using the repository browser.