Geant4  10.01.p02
G4VMscModel.hh
Go to the documentation of this file.
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: G4VMscModel.hh 80189 2014-04-07 10:16:26Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4VMscModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 07.03.2008
38 //
39 // Modifications:
40 // 07.04.2009 V.Ivanchenko moved msc methods from G4VEmModel to G4VMscModel
41 // 26.03.2012 V.Ivanchenko added transport x-section pointer and Get?Set methods
42 //
43 // Class Description:
44 //
45 // General interface to msc models
46 
47 // -------------------------------------------------------------------
48 //
49 #ifndef G4VMscModel_h
50 #define G4VMscModel_h 1
51 
52 #include <CLHEP/Units/SystemOfUnits.h>
53 
54 #include "G4VEmModel.hh"
55 #include "G4MscStepLimitType.hh"
56 #include "globals.hh"
57 #include "G4ThreeVector.hh"
58 #include "G4Track.hh"
59 #include "G4SafetyHelper.hh"
60 #include "G4VEnergyLossProcess.hh"
61 #include "G4LossTableManager.hh"
62 #include "G4PhysicsTable.hh"
63 #include "G4ThreeVector.hh"
64 #include <vector>
65 
67 
68 class G4VMscModel : public G4VEmModel
69 {
70 
71 public:
72 
73  G4VMscModel(const G4String& nam);
74 
75  virtual ~G4VMscModel();
76 
77  virtual G4double ComputeTruePathLengthLimit(const G4Track& track,
78  G4double& stepLimit);
79 
80  virtual G4double ComputeGeomPathLength(G4double truePathLength);
81 
82  virtual G4double ComputeTrueStepLength(G4double geomPathLength);
83 
85  G4double safety);
86 
87  // empty method
88  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
89  const G4MaterialCutsCouple*,
90  const G4DynamicParticle*,
91  G4double tmin,
92  G4double tmax);
93 
94  //================================================================
95  // Set parameters of multiple scattering models
96  //================================================================
97 
99 
100  inline void SetLateralDisplasmentFlag(G4bool val);
101 
102  inline void SetRangeFactor(G4double);
103 
104  inline void SetGeomFactor(G4double);
105 
106  inline void SetSkin(G4double);
107 
108  inline void SetSampleZ(G4bool);
109 
110  //================================================================
111  // Get/Set access to Physics Tables
112  //================================================================
113 
114  inline G4VEnergyLossProcess* GetIonisation() const;
115 
116  inline void SetIonisation(G4VEnergyLossProcess*,
117  const G4ParticleDefinition* part);
118 
119  //================================================================
120  // Run time methods
121  //================================================================
122 
123 protected:
124 
125  // initialisation of the ParticleChange for the model
126  // initialisation of interface with geometry and ionisation
129 
130  // convert true length to geometry length
131  inline G4double ConvertTrueToGeom(G4double& tLength, G4double& gLength);
132 
133 public:
134 
135  // compute safety
137  G4double limit= DBL_MAX);
138 
139  // compute linear distance to a geometry boundary
140  inline G4double ComputeGeomLimit(const G4Track&, G4double& presafety,
141  G4double limit);
142 
143  inline G4double GetDEDX(const G4ParticleDefinition* part,
144  G4double kineticEnergy,
145  const G4MaterialCutsCouple* couple);
146 
147  inline G4double GetRange(const G4ParticleDefinition* part,
148  G4double kineticEnergy,
149  const G4MaterialCutsCouple* couple);
150 
151  inline G4double GetEnergy(const G4ParticleDefinition* part,
152  G4double range,
153  const G4MaterialCutsCouple* couple);
154 
155  // G4MaterialCutsCouple should be defined before call to this method
156  inline
158  G4double kinEnergy);
159 
160 private:
161 
162  // hide assignment operator
164  G4VMscModel(const G4VMscModel&);
165 
170 
174 
175 protected:
176 
185 
188 
191 
192 };
193 
194 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196 
198 {
199  latDisplasment = val;
200 }
201 
202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
203 
205 {
206  skin = val;
207 }
208 
209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210 
212 {
213  facrange = val;
214 }
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217 
219 {
220  facgeom = val;
221 }
222 
223 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
224 
226 {
227  steppingAlgorithm = val;
228 }
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
231 
233 {
234  samplez = val;
235 }
236 
237 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
238 
240  G4double limit)
241 {
242  return safetyHelper->ComputeSafety(position, limit);
243 }
244 
245 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
246 
248  G4double& glength)
249 {
250  glength = ComputeGeomPathLength(tlength);
251  // should return true length
252  return tlength;
253 }
254 
255 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
256 
258  G4double& presafety,
259  G4double limit)
260 {
261  G4double res = geomMax;
262  if(track.GetVolume() != safetyHelper->GetWorldVolume()) {
264  track.GetStep()->GetPreStepPoint()->GetPosition(),
265  track.GetMomentumDirection(),
266  limit, presafety);
267  }
268  return res;
269 }
270 
271 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
272 
273 inline G4double
275  G4double kinEnergy, const G4MaterialCutsCouple* couple)
276 {
277  G4double x;
278  if(ionisation) { x = ionisation->GetDEDX(kinEnergy, couple); }
279  else {
280  G4double q = part->GetPDGCharge()*inveplus;
281  x = dedx*q*q;
282  }
283  return x;
284 }
285 
286 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
287 
288 inline G4double
290  G4double kinEnergy, const G4MaterialCutsCouple* couple)
291 {
292  //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " "
293  // << ionisation << " " << part->GetParticleName()
294  // << G4endl;
295  localtkin = kinEnergy;
296  if(ionisation) {
297  localrange = ionisation->GetRangeForLoss(kinEnergy, couple);
298  } else {
299  G4double q = part->GetPDGCharge()*inveplus;
300  localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
301  }
302  //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
303  return localrange;
304 }
305 
306 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
307 
308 inline G4double
310  G4double range, const G4MaterialCutsCouple* couple)
311 {
312  G4double e;
313  //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << " " << ionisation
314  // << " Rlocal(mm)= " << localrange << " Elocal(MeV)= " << localtkin
315  // << G4endl;
316  if(ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
317  else {
318  e = localtkin;
319  if(localrange > range) {
320  G4double q = part->GetPDGCharge()*inveplus;
321  e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
322  }
323  }
324  return e;
325 }
326 
327 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
328 
330 {
331  return ionisation;
332 }
333 
334 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
335 
337  const G4ParticleDefinition* part)
338 {
339  ionisation = p;
340  currentPart = part;
341 }
342 
343 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
344 
345 inline G4double
347  G4double ekin)
348 {
349  G4double x;
350  if(xSectionTable) {
351  G4int idx = CurrentCouple()->GetIndex();
352  x = (*xSectionTable)[(*theDensityIdx)[idx]]->Value(ekin, idxTable)
353  *(*theDensityFactor)[idx]/(ekin*ekin);
354  } else {
355  x = CrossSectionPerVolume(CurrentCouple()->GetMaterial(), part, ekin,
356  0.0, DBL_MAX);
357  }
358  if(0.0 >= x) { x = DBL_MAX; }
359  else { x = 1.0/x; }
360  return x;
361 }
362 
363 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
364 
365 #endif
366 
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4double facgeom
Definition: G4VMscModel.hh:178
size_t idxTable
Definition: G4VEmModel.hh:407
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:258
CLHEP::Hep3Vector G4ThreeVector
G4double dtrl
Definition: G4VMscModel.hh:181
G4double facrange
Definition: G4VMscModel.hh:177
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4double lambdalimit
Definition: G4VMscModel.hh:182
G4bool samplez
Definition: G4VMscModel.hh:189
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
Definition: G4VMscModel.hh:247
G4double skin
Definition: G4VMscModel.hh:180
G4double dedx
Definition: G4VMscModel.hh:171
G4double GetDensity() const
Definition: G4Material.hh:180
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:197
G4bool latDisplasment
Definition: G4VMscModel.hh:190
G4double geomMax
Definition: G4VMscModel.hh:184
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:225
G4double localtkin
Definition: G4VMscModel.hh:172
const G4Step * GetStep() const
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)
Definition: G4VMscModel.hh:239
int G4int
Definition: G4Types.hh:78
G4StepPoint * GetPreStepPoint() const
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:309
G4VEnergyLossProcess * GetIonisation() const
Definition: G4VMscModel.hh:329
virtual G4double ComputeTrueStepLength(G4double geomPathLength)
Definition: G4VMscModel.cc:160
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:434
#define position
Definition: xmlparse.cc:605
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)
Definition: G4VMscModel.cc:146
void SetSampleZ(G4bool)
Definition: G4VMscModel.hh:232
const G4ThreeVector & GetPosition() const
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:289
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=0)
Definition: G4VMscModel.cc:90
bool G4bool
Definition: G4Types.hh:79
G4double geomMin
Definition: G4VMscModel.hh:183
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:211
static const G4double inveplus
Definition: G4VEmModel.hh:408
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
Definition: G4VMscModel.hh:257
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:186
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:346
G4VPhysicalVolume * GetWorldVolume()
G4double localrange
Definition: G4VMscModel.hh:173
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:336
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * currentPart
Definition: G4VMscModel.hh:168
G4double facsafety
Definition: G4VMscModel.hh:179
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax)
Definition: G4VMscModel.cc:167
G4VEnergyLossProcess * ionisation
Definition: G4VMscModel.hh:167
G4double GetDEDX(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:274
G4MscStepLimitType
G4double CheckNextStep(const G4ThreeVector &position, const G4ThreeVector &direction, const G4double currentMaxStep, G4double &newSafety)
G4VMscModel & operator=(const G4VMscModel &right)
G4VPhysicalVolume * GetVolume() const
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:187
void SetGeomFactor(G4double)
Definition: G4VMscModel.hh:218
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:404
double G4double
Definition: G4Types.hh:76
void SetSkin(G4double)
Definition: G4VMscModel.hh:204
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4VMscModel(const G4String &nam)
Definition: G4VMscModel.cc:58
G4SafetyHelper * safetyHelper
Definition: G4VMscModel.hh:166
virtual G4double ComputeGeomPathLength(G4double truePathLength)
Definition: G4VMscModel.cc:153
virtual ~G4VMscModel()
Definition: G4VMscModel.cc:84
const G4Material * GetMaterial() const
G4LossTableManager * man
Definition: G4VMscModel.hh:169
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
Definition: G4VMscModel.cc:139