source: trunk/source/processes/electromagnetic/standard/include/G4PAIModel.hh @ 1058

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

maj sur la beta de geant 4.9.3

File size: 6.3 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: G4PAIModel.hh,v 1.22 2009/02/19 19:17:50 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name:     G4PAIModel
35//
36// Author:        V. Grichine based on Vladimir Ivanchenko  code
37//
38// Creation date: 05.10.2003
39//
40// Modifications:
41// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
42// 26-09-07 Fixed tmax computation (V.Ivantchenko)
43//
44//
45// Class Description:
46//
47// Implementation of PAI model of energy loss and
48// delta-electron production by heavy charged particles
49
50// -------------------------------------------------------------------
51//
52
53#ifndef G4PAIModel_h
54#define G4PAIModel_h 1
55
56#include <vector>
57#include "G4VEmModel.hh"
58#include "globals.hh"
59#include "G4VEmFluctuationModel.hh"
60#include "G4PAIySection.hh"
61
62class G4PhysicsLogVector;
63class G4PhysicsTable;
64class G4Region;
65class G4MaterialCutsCouple;
66class G4ParticleChangeForLoss;
67
68class G4PAIModel : public G4VEmModel, public G4VEmFluctuationModel
69{
70
71public:
72
73  G4PAIModel(const G4ParticleDefinition* p = 0, const G4String& nam = "PAI");
74
75  virtual ~G4PAIModel();
76
77  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
78
79  virtual void InitialiseMe(const G4ParticleDefinition*);
80
81  virtual G4double ComputeDEDXPerVolume(const G4Material*,
82                               const G4ParticleDefinition*,
83                               G4double kineticEnergy,
84                               G4double cutEnergy);
85
86  virtual G4double CrossSectionPerVolume(const G4Material*,
87                                const G4ParticleDefinition*,
88                                G4double kineticEnergy,
89                                G4double cutEnergy,
90                                G4double maxEnergy);
91
92  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
93                                 const G4MaterialCutsCouple*,
94                                 const G4DynamicParticle*,
95                                 G4double tmin,
96                                 G4double maxEnergy);
97
98  virtual G4double SampleFluctuations(const G4Material*,
99                                      const G4DynamicParticle*,
100                                      G4double&,
101                                      G4double&,
102                                      G4double&);
103
104  virtual G4double Dispersion(    const G4Material*,
105                                  const G4DynamicParticle*,
106                                  G4double&,
107                                  G4double&);
108
109  void     DefineForRegion(const G4Region* r) ;
110  void     ComputeSandiaPhotoAbsCof();
111  void     BuildPAIonisationTable();
112  void     BuildLambdaVector();
113
114  G4double GetdNdxCut( G4int iPlace, G4double transferCut);
115  G4double GetdEdxCut( G4int iPlace, G4double transferCut);
116  G4double GetPostStepTransfer( G4double scaledTkin );
117  G4double GetEnergyTransfer( G4int iPlace,
118                              G4double position,
119                              G4int iTransfer );
120
121  void SetVerboseLevel(G4int verbose){fVerbose=verbose;};
122
123protected:
124
125  G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
126                                    G4double kinEnergy);
127
128private:
129
130  void SetParticle(const G4ParticleDefinition* p);
131
132  // hide assignment operator
133  G4PAIModel & operator=(const  G4PAIModel &right);
134  G4PAIModel(const  G4PAIModel&);
135
136  // The vector over proton kinetic energies: the range of gammas
137  G4int                fVerbose; 
138  G4double             fLowestGamma;
139  G4double             fHighestGamma;
140  G4double             fLowestKineticEnergy;
141  G4double             fHighestKineticEnergy;
142  G4int                fTotBin;
143  G4int                fMeanNumber;
144  G4PhysicsLogVector*  fParticleEnergyVector ;
145  G4PAIySection        fPAIySection;
146
147  // vectors
148
149  G4PhysicsTable*                    fPAItransferTable;
150  std::vector<G4PhysicsTable*>       fPAIxscBank;
151
152  G4PhysicsTable*                    fPAIdEdxTable;
153  std::vector<G4PhysicsTable*>       fPAIdEdxBank;
154
155  std::vector<const G4MaterialCutsCouple*> fMaterialCutsCoupleVector;
156  std::vector<const G4Region*>       fPAIRegionVector;
157
158  const G4MaterialCutsCouple*        fCutCouple;
159  const G4Material*                  fMaterial;
160  G4double                           fDeltaCutInKinEnergy; 
161
162  size_t                             fMatIndex ; 
163  G4double**                         fSandiaPhotoAbsCof ;
164  G4int                              fSandiaIntervalNumber ;
165
166  G4PhysicsLogVector*                fdEdxVector ;
167  std::vector<G4PhysicsLogVector*>   fdEdxTable ;
168
169  G4PhysicsLogVector*                fLambdaVector ;
170  std::vector<G4PhysicsLogVector*>   fLambdaTable ;
171
172  G4PhysicsLogVector*                fdNdxCutVector ;
173  std::vector<G4PhysicsLogVector*>   fdNdxCutTable ;
174
175  const G4ParticleDefinition* fParticle;
176  const G4ParticleDefinition* fElectron;
177  const G4ParticleDefinition* fPositron;
178  G4ParticleChangeForLoss*    fParticleChange;
179
180  G4double fMass;
181  G4double fSpin;
182  G4double fChargeSquare;
183  G4double fRatio;
184  G4double fHighKinEnergy;
185  G4double fLowKinEnergy;
186  G4double fTwoln10;
187  G4double fBg2lim; 
188  G4double fTaulim;
189  G4double fQc;
190
191  G4bool   isInitialised;
192};
193
194#endif
195
196
197
198
199
200
201
Note: See TracBrowser for help on using the repository browser.