Geant4_10
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 71484 2013-06-17 08:12:51Z 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 
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
136  inline G4double ComputeSafety(const G4ThreeVector& position, G4double limit);
137 
138  // compute linear distance to a geometry boundary
139  inline G4double ComputeGeomLimit(const G4Track&, G4double& presafety,
140  G4double limit);
141 
142  inline G4double GetDEDX(const G4ParticleDefinition* part,
143  G4double kineticEnergy,
144  const G4MaterialCutsCouple* couple);
145 
146  inline G4double GetRange(const G4ParticleDefinition* part,
147  G4double kineticEnergy,
148  const G4MaterialCutsCouple* couple);
149 
150  inline G4double GetEnergy(const G4ParticleDefinition* part,
151  G4double range,
152  const G4MaterialCutsCouple* couple);
153 
154  // G4MaterialCutsCouple should be defined before call to this method
155  inline
157  G4double kinEnergy);
158 
159 private:
160 
161  // hide assignment operator
162  G4VMscModel & operator=(const G4VMscModel &right);
163  G4VMscModel(const G4VMscModel&);
164 
165  G4SafetyHelper* safetyHelper;
166  G4VEnergyLossProcess* ionisation;
167  const G4ParticleDefinition* currentPart;
168  G4LossTableManager* man;
169 
170  G4double dedx;
171  G4double localtkin;
172  G4double localrange;
173 
174 protected:
175 
184 
187 
190 
191 };
192 
193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
194 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195 
197 {
198  latDisplasment = val;
199 }
200 
201 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202 
204 {
205  skin = val;
206 }
207 
208 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
209 
211 {
212  facrange = val;
213 }
214 
215 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
216 
218 {
219  facgeom = val;
220 }
221 
222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
223 
225 {
226  steppingAlgorithm = val;
227 }
228 
229 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
230 
232 {
233  samplez = val;
234 }
235 
236 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237 
239  G4double)
240 {
241  return safetyHelper->ComputeSafety(position);
242 }
243 
244 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
245 
247  G4double& glength)
248 {
249  glength = ComputeGeomPathLength(tlength);
250  // should return true length
251  return tlength;
252 }
253 
254 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
255 
257  G4double& presafety,
258  G4double limit)
259 {
260  G4double res = geomMax;
261  if(track.GetVolume() != safetyHelper->GetWorldVolume()) {
262  res = safetyHelper->CheckNextStep(
263  track.GetStep()->GetPreStepPoint()->GetPosition(),
264  track.GetMomentumDirection(),
265  limit, presafety);
266  }
267  return res;
268 }
269 
270 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
271 
272 inline G4double
274  G4double kinEnergy, const G4MaterialCutsCouple* couple)
275 {
276  G4double x;
277  if(ionisation) { x = ionisation->GetDEDX(kinEnergy, couple); }
278  else {
279  G4double q = part->GetPDGCharge()/CLHEP::eplus;
280  x = dedx*q*q;
281  }
282  return x;
283 }
284 
285 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
286 
287 inline G4double
289  G4double kinEnergy, const G4MaterialCutsCouple* couple)
290 {
291  //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " " << ionisation
292  // << " " << part->GetParticleName()
293  // << G4endl;
294  localtkin = kinEnergy;
295  if(ionisation) {
296  localrange = ionisation->GetRangeForLoss(kinEnergy, couple);
297  } else {
298  G4double q = part->GetPDGCharge()/CLHEP::eplus;
299  localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
300  }
301  //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
302  return localrange;
303 }
304 
305 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
306 
307 inline G4double
309  G4double range, const G4MaterialCutsCouple* couple)
310 {
311  G4double e;
312  //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << " " << ionisation
313  // << " Rlocal(mm)= " << localrange << " Elocal(MeV)= " << localtkin
314  // << G4endl;
315  if(ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
316  else {
317  e = localtkin;
318  if(localrange > range) {
319  G4double q = part->GetPDGCharge()/CLHEP::eplus;
320  e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
321  }
322  }
323  return e;
324 }
325 
326 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
327 
329 {
330  return ionisation;
331 }
332 
333 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
334 
336  const G4ParticleDefinition* part)
337 {
338  ionisation = p;
339  currentPart = part;
340 }
341 
342 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
343 
344 inline G4double
346  G4double ekin)
347 {
348  G4double x;
349  if(xSectionTable) {
350  G4int idx = CurrentCouple()->GetIndex();
351  x = (*xSectionTable)[(*theDensityIdx)[idx]]->Value(ekin, idxTable)
352  *(*theDensityFactor)[idx]/(ekin*ekin);
353  } else {
354  x = CrossSectionPerVolume(CurrentCouple()->GetMaterial(), part, ekin,
355  0.0, DBL_MAX);
356  }
357  if(0.0 >= x) { x = DBL_MAX; }
358  else { x = 1.0/x; }
359  return x;
360 }
361 
362 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
363 
364 #endif
365 
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4double facgeom
Definition: G4VMscModel.hh:177
size_t idxTable
Definition: G4VEmModel.hh:402
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
G4double dtrl
Definition: G4VMscModel.hh:180
G4double facrange
Definition: G4VMscModel.hh:176
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
const char * p
Definition: xmltok.h:285
G4double lambdalimit
Definition: G4VMscModel.hh:181
G4bool samplez
Definition: G4VMscModel.hh:188
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
Definition: G4VMscModel.hh:246
G4double skin
Definition: G4VMscModel.hh:179
G4double GetDensity() const
Definition: G4Material.hh:178
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:196
G4bool latDisplasment
Definition: G4VMscModel.hh:189
G4double geomMax
Definition: G4VMscModel.hh:183
tuple x
Definition: test.py:50
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:224
const G4Step * GetStep() const
int G4int
Definition: G4Types.hh:78
G4StepPoint * GetPreStepPoint() const
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:308
G4VEnergyLossProcess * GetIonisation() const
Definition: G4VMscModel.hh:328
virtual G4double ComputeTrueStepLength(G4double geomPathLength)
Definition: G4VMscModel.cc:159
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)
Definition: G4VMscModel.cc:145
void SetSampleZ(G4bool)
Definition: G4VMscModel.hh:231
const G4ParticleDefinition const G4Material *G4double range
const G4ThreeVector & GetPosition() const
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:288
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=0)
Definition: G4VMscModel.cc:89
bool G4bool
Definition: G4Types.hh:79
G4double geomMin
Definition: G4VMscModel.hh:182
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:210
TString part[npart]
Definition: Style.C:32
G4double ComputeGeomLimit(const G4Track &, G4double &presafety, G4double limit)
Definition: G4VMscModel.hh:256
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:185
G4double GetTransportMeanFreePath(const G4ParticleDefinition *part, G4double kinEnergy)
Definition: G4VMscModel.hh:345
G4VPhysicalVolume * GetWorldVolume()
G4double ComputeSafety(const G4ThreeVector &position, G4double limit)
Definition: G4VMscModel.hh:238
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:335
const G4ThreeVector & GetMomentumDirection() const
G4double facsafety
Definition: G4VMscModel.hh:178
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax)
Definition: G4VMscModel.cc:166
G4double GetDEDX(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:273
G4MscStepLimitType
G4double CheckNextStep(const G4ThreeVector &position, const G4ThreeVector &direction, const G4double currentMaxStep, G4double &newSafety)
G4VPhysicalVolume * GetVolume() const
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:186
void SetGeomFactor(G4double)
Definition: G4VMscModel.hh:217
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:399
double G4double
Definition: G4Types.hh:76
void SetSkin(G4double)
Definition: G4VMscModel.hh:203
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4VMscModel(const G4String &nam)
Definition: G4VMscModel.cc:58
virtual G4double ComputeGeomPathLength(G4double truePathLength)
Definition: G4VMscModel.cc:152
virtual ~G4VMscModel()
Definition: G4VMscModel.cc:83
const G4Material * GetMaterial() const
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
Definition: G4VMscModel.cc:138