source: trunk/source/processes/electromagnetic/standard/include/G4WentzelVIModel.hh

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

tag geant4.9.4 beta 1 + modifs locales

File size: 6.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: G4WentzelVIModel.hh,v 1.29 2010/05/27 14:22:05 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-01 $
28//
29// -------------------------------------------------------------------
30//
31//
32// GEANT4 Class header file
33//
34//
35// File name:     G4WentzelVIModel
36//
37// Author:        V.Ivanchenko
38//
39// Creation date: 09.04.2008 from G4MuMscModel
40//
41// Modifications:
42// 27-05-2010 V.Ivanchenko added G4WentzelOKandVIxSection class to
43//              compute cross sections and sample scattering angle
44//
45// Class Description:
46//
47// Implementation of the model of multiple scattering based on
48// G.Wentzel, Z. Phys. 40 (1927) 590.
49// H.W.Lewis, Phys Rev 78 (1950) 526.
50// J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
51// L.Urban, CERN-OPEN-2006-077.
52
53// -------------------------------------------------------------------
54//
55
56#ifndef G4WentzelVIModel_h
57#define G4WentzelVIModel_h 1
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
61#include "G4VMscModel.hh"
62#include "G4PhysicsTable.hh"
63#include "G4MaterialCutsCouple.hh"
64#include "G4WentzelOKandVIxSection.hh"
65
66class G4ParticleDefinition;
67class G4LossTableManager;
68class G4Pow;
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71
72class G4WentzelVIModel : public G4VMscModel
73{
74
75public:
76
77  G4WentzelVIModel(const G4String& nam = "WentzelVIUni");
78
79  virtual ~G4WentzelVIModel();
80
81  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
82
83  virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
84                                              G4double KineticEnergy,
85                                              G4double AtomicNumber,
86                                              G4double AtomicWeight=0., 
87                                              G4double cut = DBL_MAX,
88                                              G4double emax= DBL_MAX);
89
90  virtual void SampleScattering(const G4DynamicParticle*, G4double safety);
91
92  virtual G4double ComputeTruePathLengthLimit(const G4Track& track,
93                                              G4PhysicsTable* theLambdaTable,
94                                              G4double currentMinimalStep);
95
96  virtual G4double ComputeGeomPathLength(G4double truePathLength);
97
98  virtual G4double ComputeTrueStepLength(G4double geomStepLength);
99
100private:
101
102  G4double ComputeXSectionPerVolume();
103
104  inline G4double GetLambda(G4double kinEnergy);
105
106  inline void SetupParticle(const G4ParticleDefinition*);
107
108  inline void DefineMaterial(const G4MaterialCutsCouple*);
109
110  //  hide assignment operator
111  G4WentzelVIModel & operator=(const  G4WentzelVIModel &right);
112  G4WentzelVIModel(const  G4WentzelVIModel&);
113
114  G4LossTableManager*       theManager;
115  G4ParticleChangeForMSC*   fParticleChange;
116  G4WentzelOKandVIxSection* wokvi;
117  G4Pow*                    fG4pow;
118
119  G4PhysicsTable*           theLambdaTable;
120  const G4DataVector*       currentCuts;
121
122  G4double tlimitminfix;
123  G4double invsqrt12;
124
125  // cache kinematics
126  G4double preKinEnergy;
127  G4double tPathLength;
128  G4double zPathLength;
129  G4double lambdaeff;
130  G4double currentRange; 
131
132  // data for single scattering mode
133  G4double xtsec;
134  std::vector<G4double> xsecn;
135  std::vector<G4double> prob;
136  G4int    nelments;
137
138  G4double numlimit;
139
140  // cache material
141  G4int    currentMaterialIndex;
142  const G4MaterialCutsCouple* currentCouple;
143  const G4Material* currentMaterial;
144
145  // single scattering parameters
146  G4double cosThetaMin;
147  G4double cosThetaMax;
148  G4double cosTetMaxNuc;
149
150  // projectile
151  const G4ParticleDefinition* particle;
152  G4double lowEnergyLimit;
153
154  // flags
155  G4bool   isInitialized;
156  G4bool   inside;
157};
158
159//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161
162inline
163void G4WentzelVIModel::DefineMaterial(const G4MaterialCutsCouple* cup) 
164{ 
165  if(cup != currentCouple) {
166    currentCouple = cup;
167    currentMaterial = cup->GetMaterial();
168    currentMaterialIndex = currentCouple->GetIndex(); 
169  }
170}
171
172//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
173
174inline
175G4double G4WentzelVIModel::GetLambda(G4double e)
176{
177  G4double x;
178  if(theLambdaTable) {
179    x = ((*theLambdaTable)[currentMaterialIndex])->Value(e);
180  } else {
181    x = CrossSection(currentCouple,particle,e,
182                     (*currentCuts)[currentMaterialIndex]);
183  }
184  if(x > DBL_MIN) { x = 1./x; }
185  else            { x = DBL_MAX; }
186  return x;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
190
191inline void G4WentzelVIModel::SetupParticle(const G4ParticleDefinition* p)
192{
193  // Initialise mass and charge
194  if(p != particle) {
195    particle = p;
196    wokvi->SetupParticle(p);
197  }
198}
199
200//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
201
202#endif
203
Note: See TracBrowser for help on using the repository browser.