source: trunk/source/processes/electromagnetic/utils/include/G4EmCalculator.hh @ 961

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

update processes

File size: 10.0 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: G4EmCalculator.hh,v 1.18 2007/03/15 12:34:46 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29//
30// -------------------------------------------------------------------
31//
32// GEANT4 Class header file
33//
34//
35// File name:     G4EmCalculator
36//
37// Author:        Vladimir Ivanchenko
38//
39// Creation date: 27.06.2004
40//
41// Modifications:
42// 17.11.2004 Change signature of methods, add new methods (V.Ivanchenko)
43// 11.01.2006 Add GetCSDARange (V.Ivanchenko)
44// 26.01.2006 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
45// 22.03.2006 Add ComputeElectronicDEDX and ComputeTotalDEDX (V.Ivanchenko)
46// 29.09.2006 Add member loweModel (V.Ivanchenko)
47// 15.03.2007 Add ComputeEnergyCutFromRangeCut methods (V.Ivanchenko)
48//
49// Class Description:
50//
51// Provide access to dE/dx and cross sections
52
53// -------------------------------------------------------------------
54//
55
56#ifndef G4EmCalculator_h
57#define G4EmCalculator_h 1
58
59#include <vector>
60#include "globals.hh"
61#include "G4DataVector.hh"
62
63class G4LossTableManager;
64class G4Material;
65class G4MaterialCutsCouple;
66class G4ParticleDefinition;
67class G4PhysicsTable;
68class G4VEmModel;
69class G4VEnergyLossProcess;
70class G4ionEffectiveCharge;
71class G4Region;
72class G4Element;
73class G4EmCorrections;
74
75class G4EmCalculator
76{
77
78public:
79
80  G4EmCalculator();
81
82  ~G4EmCalculator();
83
84  // Methods to access precalculated dE/dx and cross sections
85  // Materials should exist in the list of the G4MaterialCutsCouple
86
87  G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition*, const G4Material*,
88                   const G4Region* r = 0);
89  G4double GetDEDX(G4double kinEnergy, const G4String& part, const G4String& mat,
90                   const G4String& s = "world");
91
92  G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition*, 
93                                     const G4Material*,
94                                     const G4Region* r = 0);
95  G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4String& part, 
96                                     const G4String& mat,
97                                     const G4String& s = "world");
98
99  G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition*, 
100                        const G4Material*,
101                        const G4Region* r = 0);
102  G4double GetCSDARange(G4double kinEnergy, const G4String& part, 
103                        const G4String& mat,
104                        const G4String& s = "world");
105
106  G4double GetRange(G4double kinEnergy, const G4ParticleDefinition*, 
107                        const G4Material*,
108                        const G4Region* r = 0);
109  G4double GetRange(G4double kinEnergy, const G4String& part, 
110                        const G4String& mat,
111                        const G4String& s = "world");
112
113  G4double GetKinEnergy(G4double range, const G4ParticleDefinition*, const G4Material*,
114                   const G4Region* r = 0);
115  G4double GetKinEnergy(G4double range, const G4String& part, const G4String& mat,
116                   const G4String& s = "world");
117
118  G4double GetCrossSectionPerVolume(
119                   G4double kinEnergy, const G4ParticleDefinition*,
120                   const G4String& processName,  const G4Material*,
121                   const G4Region* r = 0);
122  G4double GetCrossSectionPerVolume(
123                   G4double kinEnergy, const G4String& part, const G4String& proc,
124                   const G4String& mat, const G4String& s = "world");
125
126  G4double GetMeanFreePath(G4double kinEnergy, const G4ParticleDefinition*,
127                   const G4String& processName,  const G4Material*,
128                   const G4Region* r = 0);
129  G4double GetMeanFreePath(G4double kinEnergy, const G4String& part, const G4String& proc,
130                   const G4String& mat, const G4String& s = "world");
131
132  void PrintDEDXTable(const G4ParticleDefinition*);
133
134  void PrintRangeTable(const G4ParticleDefinition*);
135
136  void PrintInverseRangeTable(const G4ParticleDefinition*);
137
138  // Methods to calculate dE/dx and cross sections "on fly"
139  // Existing tables and G4MaterialCutsCouples are not used
140
141  G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition*,
142                       const G4String& processName,  const G4Material*,
143                       G4double cut = DBL_MAX);
144  G4double ComputeDEDX(G4double kinEnergy, const G4String& part, const G4String& proc,
145                       const G4String& mat, G4double cut = DBL_MAX);
146
147  G4double ComputeElectronicDEDX(G4double kinEnergy, const G4ParticleDefinition*,
148                                 const G4Material* mat, G4double cut = DBL_MAX);
149  G4double ComputeElectronicDEDX(G4double kinEnergy, const G4String& part,
150                                 const G4String& mat, G4double cut = DBL_MAX);
151
152  G4double ComputeNuclearDEDX(G4double kinEnergy, const G4ParticleDefinition*, const G4Material*);
153  G4double ComputeNuclearDEDX(G4double kinEnergy, const G4String& part, const G4String& mat);
154
155  G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition*, 
156                            const G4Material*, G4double cut = DBL_MAX);
157  G4double ComputeTotalDEDX(G4double kinEnergy, const G4String& part, 
158                            const G4String& mat, G4double cut = DBL_MAX);
159
160  G4double ComputeCrossSectionPerVolume(
161                       G4double kinEnergy, const G4ParticleDefinition*,
162                       const G4String& processName,  const G4Material*,
163                       G4double cut = 0.0);
164  G4double ComputeCrossSectionPerVolume(
165                       G4double kinEnergy, const G4String& part, const G4String& proc,
166                       const G4String& mat, G4double cut = 0.0);
167
168  G4double ComputeCrossSectionPerAtom(
169                       G4double kinEnergy, const G4ParticleDefinition*,
170                       const G4String& processName, G4double Z, G4double A,
171                       G4double cut = 0.0);
172  G4double ComputeCrossSectionPerAtom(
173                       G4double kinEnergy, const G4String& part,
174                       const G4String& processName, const G4Element*,
175                       G4double cut = 0.0);
176
177  G4double ComputeMeanFreePath(
178                       G4double kinEnergy, const G4ParticleDefinition*,
179                       const G4String& processName,  const G4Material*,
180                       G4double cut = 0.0);
181  G4double ComputeMeanFreePath(
182                       G4double kinEnergy, const G4String&, const G4String&,
183                       const G4String& processName, G4double cut = 0.0);
184
185  G4double ComputeEnergyCutFromRangeCut(
186                       G4double range, const G4ParticleDefinition*,
187                       const G4Material*);
188  G4double ComputeEnergyCutFromRangeCut(
189                       G4double range, const G4String&,
190                       const G4String&);
191
192  const G4ParticleDefinition* FindParticle(const G4String&);
193
194  const G4Material* FindMaterial(const G4String&);
195
196  const G4Region* FindRegion(const G4String&);
197
198  const G4MaterialCutsCouple* FindCouple(const G4Material*, const G4Region* r = 0);
199
200  void SetVerbose(G4int val);
201
202private:
203
204  G4bool UpdateParticle(const G4ParticleDefinition*, G4double kinEnergy);
205
206  G4bool UpdateCouple(const G4Material*, G4double cut);
207
208  void FindLambdaTable(const G4ParticleDefinition*, const G4String& processName);
209
210  G4bool FindEmModel(const G4ParticleDefinition*, 
211                     const G4String& processName,
212                           G4double kinEnergy);
213
214  G4VEnergyLossProcess* FindEnergyLossProcess(const G4ParticleDefinition*);
215
216  G4EmCalculator & operator=(const  G4EmCalculator &right);
217  G4EmCalculator(const  G4EmCalculator&);
218
219  std::vector<const G4Material*>            localMaterials;
220  std::vector<const G4MaterialCutsCouple*>  localCouples;
221
222  G4LossTableManager*          manager;
223  G4EmCorrections*             corr; 
224  G4DataVector                 localCuts;
225  G4int                        nLocalMaterials;
226
227  G4int                        verbose;
228
229  // cash
230  G4int                        currentCoupleIndex;
231  const G4MaterialCutsCouple*  currentCouple;
232  const G4Material*            currentMaterial;
233  const G4ParticleDefinition*  currentParticle;
234  const G4ParticleDefinition*  baseParticle;
235  const G4PhysicsTable*        currentLambda;
236        G4VEmModel*            currentModel;
237        G4VEmModel*            loweModel;
238        G4VEnergyLossProcess*  currentProcess;
239
240  const G4ParticleDefinition*  theGenericIon;
241  G4ionEffectiveCharge*        ionEffCharge;
242
243  G4String                     currentName;
244  G4double                     currentCut;
245  G4double                     chargeSquare;
246  G4double                     massRatio;
247  G4double                     mass;
248  G4bool                       isIon;
249  G4bool                       isApplicable;
250
251  G4String                     currentParticleName;
252  G4String                     currentMaterialName;
253  G4String                     currentProcessName;
254};
255
256//....oooOO0OOooo.......oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
257
258#endif
Note: See TracBrowser for help on using the repository browser.