source: trunk/source/processes/electromagnetic/lowenergy/include/G4PenelopeBremsstrahlungModel.hh @ 1350

Last change on this file since 1350 was 1347, checked in by garnier, 14 years ago

geant4 tag 9.4

File size: 5.1 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: G4PenelopeBremsstrahlungModel.hh,v 1.2 2009/04/17 10:29:20 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-ref-00 $
28//
29// Author: Luciano Pandola
30//
31// History:
32// -----------
33// 05 Dec 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- bremsstrahlung
40// with Penelope Model
41// -------------------------------------------------------------------
42
43#ifndef G4PENELOPEBREMSSTRAHLUNGMODEL_HH
44#define G4PENELOPEBREMSSTRAHLUNGMODEL_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 G4VEnergySpectrum;
59class G4PenelopeBremsstrahlungAngular;
60class G4PenelopeBremsstrahlungContinuous;
61
62class G4PenelopeBremsstrahlungModel : public G4VEmModel
63{
64
65public:
66 
67  G4PenelopeBremsstrahlungModel(const G4ParticleDefinition* p=0,
68                         const G4String& processName ="PenelopeBrem");
69 
70  virtual ~G4PenelopeBremsstrahlungModel();
71
72  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
73
74  virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
75                                              G4double kinEnergy,
76                                              G4double Z,
77                                              G4double A=0,
78                                              G4double cut=0,
79                                              G4double emax=DBL_MAX);
80                                         
81  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
82                                 const G4MaterialCutsCouple*,
83                                 const G4DynamicParticle*,
84                                 G4double tmin,
85                                 G4double maxEnergy);
86                                   
87  virtual G4double ComputeDEDXPerVolume(const G4Material*,
88                               const G4ParticleDefinition*,
89                               G4double kineticEnergy,
90                               G4double cutEnergy);
91                                 
92  // min cut in kinetic energy allowed by the model
93  virtual G4double MinEnergyCut(const G4ParticleDefinition*,
94                                const G4MaterialCutsCouple*);
95
96  void SetVerbosityLevel(G4int lev){verboseLevel = lev;};
97  G4int GetVerbosityLevel(){return verboseLevel;};
98
99
100protected:
101  G4ParticleChangeForLoss* fParticleChange;
102
103private:
104 
105  G4PenelopeBremsstrahlungModel & operator=(const G4PenelopeBremsstrahlungModel &right);
106  G4PenelopeBremsstrahlungModel(const G4PenelopeBremsstrahlungModel&);
107
108 
109  //Intrinsic energy limits of the model: cannot be extended by the parent
110  // process
111  G4double fIntrinsicLowEnergyLimit;
112  G4double fIntrinsicHighEnergyLimit;
113
114  G4int verboseLevel;
115
116  G4bool isInitialised;
117
118  G4VEnergySpectrum* energySpectrum;
119 
120
121  G4PenelopeBremsstrahlungAngular* GetAngularDataForZ(G4int iZ);
122  // Map to the objects containing tha angular data
123  std::map<G4int,G4PenelopeBremsstrahlungAngular*> *angularData;
124
125  G4PenelopeBremsstrahlungContinuous* GetStoppingPowerData(G4int iZ,G4double energyCut,
126                                                           const G4ParticleDefinition*);
127  std::map<std::pair<G4int,G4double>,G4PenelopeBremsstrahlungContinuous*> *stoppingPowerData;
128
129
130  G4int SampleRandomAtom(const G4MaterialCutsCouple*,G4double energy) const;
131  G4VCrossSectionHandler* crossSectionHandler;
132
133 
134};
135
136#endif
137
Note: See TracBrowser for help on using the repository browser.