ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VMscModel.hh
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4VMscModel.hh
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 // -------------------------------------------------------------------
27 //
28 // GEANT4 Class header file
29 //
30 //
31 // File name: G4VMscModel
32 //
33 // Author: Vladimir Ivanchenko
34 //
35 // Creation date: 07.03.2008
36 //
37 // Modifications:
38 // 07.04.2009 V.Ivanchenko moved msc methods from G4VEmModel to G4VMscModel
39 // 26.03.2012 V.Ivanchenko added transport x-section pointer and Get?Set methods
40 //
41 // Class Description:
42 //
43 // General interface to msc models
44 
45 // -------------------------------------------------------------------
46 //
47 #ifndef G4VMscModel_h
48 #define G4VMscModel_h 1
49 
51 
52 #include "G4VEmModel.hh"
53 #include "G4MscStepLimitType.hh"
54 #include "globals.hh"
55 #include "G4ThreeVector.hh"
56 #include "G4Track.hh"
57 #include "G4SafetyHelper.hh"
58 #include "G4VEnergyLossProcess.hh"
59 #include "G4PhysicsTable.hh"
60 #include "G4ThreeVector.hh"
61 #include <vector>
62 
65 
66 class G4VMscModel : public G4VEmModel
67 {
68 
69 public:
70 
71  explicit G4VMscModel(const G4String& nam);
72 
73  ~G4VMscModel() override;
74 
76  G4double& stepLimit) = 0;
77 
78  virtual G4double ComputeGeomPathLength(G4double truePathLength) = 0;
79 
80  virtual G4double ComputeTrueStepLength(G4double geomPathLength) = 0;
81 
83  G4double safety) = 0;
84 
86 
87  // empty method
88  void SampleSecondaries(std::vector<G4DynamicParticle*>*,
89  const G4MaterialCutsCouple*,
90  const G4DynamicParticle*,
91  G4double tmin, G4double tmax) override;
92 
93  //================================================================
94  // Set parameters of multiple scattering models
95  //================================================================
96 
98 
99  inline void SetLateralDisplasmentFlag(G4bool val);
100 
101  inline void SetRangeFactor(G4double);
102 
103  inline void SetGeomFactor(G4double);
104 
105  inline void SetSkin(G4double);
106 
107  inline void SetLambdaLimit(G4double);
108 
109  inline void SetSafetyFactor(G4double);
110 
111  inline void SetSampleZ(G4bool);
112 
113  //================================================================
114  // Get/Set access to Physics Tables
115  //================================================================
116 
117  inline G4VEnergyLossProcess* GetIonisation() const;
118 
119  inline void SetIonisation(G4VEnergyLossProcess*,
120  const G4ParticleDefinition* part);
121 
122  //================================================================
123  // Run time methods
124  //================================================================
125 
126 protected:
127 
128  // initialisation of the ParticleChange for the model
129  // initialisation of interface with geometry and ionisation
132 
133  // convert true length to geometry length
134  inline G4double ConvertTrueToGeom(G4double& tLength, G4double& gLength);
135 
136 public:
137 
138  // compute safety
140  G4double limit= DBL_MAX);
141 
142  // compute linear distance to a geometry boundary
143  inline G4double ComputeGeomLimit(const G4Track&, G4double& presafety,
144  G4double limit);
145 
146  inline G4double GetDEDX(const G4ParticleDefinition* part,
147  G4double kineticEnergy,
148  const G4MaterialCutsCouple* couple);
149  inline G4double GetDEDX(const G4ParticleDefinition* part,
150  G4double kineticEnergy,
151  const G4MaterialCutsCouple* couple,
152  G4double logKineticEnergy);
153 
154  inline G4double GetRange(const G4ParticleDefinition* part,
155  G4double kineticEnergy,
156  const G4MaterialCutsCouple* couple);
157  inline G4double GetRange(const G4ParticleDefinition* part,
158  G4double kineticEnergy,
159  const G4MaterialCutsCouple* couple,
160  G4double logKineticEnergy);
161 
162  inline G4double GetEnergy(const G4ParticleDefinition* part,
163  G4double range,
164  const G4MaterialCutsCouple* couple);
165 
166  // G4MaterialCutsCouple should be defined before call to this method
167  inline
169  G4double kinEnergy);
170  inline
172  G4double kinEnergy,
173  G4double logKinEnergy);
174 
175 private:
176 
177  // hide assignment operator
178  G4VMscModel & operator=(const G4VMscModel &right) = delete;
179  G4VMscModel(const G4VMscModel&) = delete;
180 
184 
188 
189 protected:
190 
199 
202 
205 
206 };
207 
208 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210 
212 {
213  if(!IsLocked()) { latDisplasment = val; }
214 }
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217 
219 {
220  if(!IsLocked()) { skin = val; }
221 }
222 
223 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
224 
226 {
227  if(!IsLocked()) { facrange = val; }
228 }
229 
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
231 
233 {
234  if(!IsLocked()) { facgeom = val; }
235 }
236 
237 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
238 
240 {
241  if(!IsLocked()) { lambdalimit = val; }
242 }
243 
244 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
245 
247 {
248  if(!IsLocked()) { facsafety = val; }
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
252 
254 {
255  if(!IsLocked()) { steppingAlgorithm = val; }
256 }
257 
258 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
259 
261 {
262  if(!IsLocked()) { samplez = val; }
263 }
264 
265 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
266 
268  G4double limit)
269 {
270  return safetyHelper->ComputeSafety(position, limit);
271 }
272 
273 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
274 
276  G4double& glength)
277 {
278  glength = ComputeGeomPathLength(tlength);
279  // should return true length
280  return tlength;
281 }
282 
283 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
284 
286  G4double& presafety,
287  G4double limit)
288 {
289  return safetyHelper->CheckNextStep(
290  track.GetStep()->GetPreStepPoint()->GetPosition(),
291  track.GetMomentumDirection(),
292  limit, presafety);
293 }
294 
295 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
296 
297 inline G4double
299  const G4MaterialCutsCouple* couple)
300 {
301  G4double x;
302  if (ionisation) {
303  x = ionisation->GetDEDX(kinEnergy, couple);
304  } else {
305  const G4double q = part->GetPDGCharge()*inveplus;
306  x = dedx*q*q;
307  }
308  return x;
309 }
310 
311 inline G4double
313  const G4MaterialCutsCouple* couple, G4double logKinEnergy)
314 {
315  G4double x;
316  if (ionisation) {
317  x = ionisation->GetDEDX(kinEnergy, couple, logKinEnergy);
318  } else {
319  const G4double q = part->GetPDGCharge()*inveplus;
320  x = dedx*q*q;
321  }
322  return x;
323 }
324 
325 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
326 
327 inline G4double
329  const G4MaterialCutsCouple* couple)
330 {
331  //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " "
332  // << ionisation << " " << part->GetParticleName()
333  // << G4endl;
334  localtkin = kinEnergy;
335  if (ionisation) {
336  localrange = ionisation->GetRangeForLoss(kinEnergy, couple);
337  } else {
338  const G4double q = part->GetPDGCharge()*inveplus;
339  localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
340  }
341  //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
342  return localrange;
343 }
344 
345 inline G4double
347  const G4MaterialCutsCouple* couple, G4double logKinEnergy)
348 {
349  //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " "
350  // << ionisation << " " << part->GetParticleName()
351  // << G4endl;
352  localtkin = kinEnergy;
353  if (ionisation) {
354  localrange = ionisation->GetRangeForLoss(kinEnergy, couple, logKinEnergy);
355  } else {
356  const G4double q = part->GetPDGCharge()*inveplus;
357  localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
358  }
359  //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
360  return localrange;
361 }
362 
363 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
364 
365 inline G4double
367  G4double range, const G4MaterialCutsCouple* couple)
368 {
369  G4double e;
370  //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << " " << ionisation
371  // << " Rlocal(mm)= " << localrange << " Elocal(MeV)= " << localtkin
372  // << G4endl;
373  if(ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
374  else {
375  e = localtkin;
376  if(localrange > range) {
377  G4double q = part->GetPDGCharge()*inveplus;
378  e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
379  }
380  }
381  return e;
382 }
383 
384 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
385 
387 {
388  return ionisation;
389 }
390 
391 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
392 
394  const G4ParticleDefinition* part)
395 {
396  ionisation = p;
397  currentPart = part;
398 }
399 
400 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
401 
402 inline G4double
404  G4double ekin)
405 {
406  G4double x;
407  if (xSectionTable) {
408  const G4int idx = CurrentCouple()->GetIndex();
409  x = (*xSectionTable)[idx]->Value(ekin, idxTable)/(ekin*ekin);
410  } else {
411  x = CrossSectionPerVolume(CurrentCouple()->GetMaterial(), part, ekin,
412  0.0, DBL_MAX);
413  }
414  return (x > 0.0) ? 1.0/x : DBL_MAX;
415 }
416 
417 inline G4double
419  G4double ekin, G4double logekin)
420 {
421  G4double x;
422  if (xSectionTable) {
423  const G4int idx = CurrentCouple()->GetIndex();
424  x = (*xSectionTable)[idx]->LogVectorValue(ekin, logekin)/(ekin*ekin);
425  } else {
426  x = CrossSectionPerVolume(CurrentCouple()->GetMaterial(), part, ekin,
427  0.0, DBL_MAX);
428  }
429  return (x > 0.0) ? 1.0/x : DBL_MAX;
430 }
431 
432 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
433 
434 #endif
435