Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/utils/include/G4VMscModel.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ]

  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 
 50 #include <CLHEP/Units/SystemOfUnits.h>
 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 
 63 class G4ParticleChangeForMSC;
 64 class G4ParticleDefinition;
 65 
 66 class G4VMscModel : public G4VEmModel
 67 {
 68 
 69 public:
 70 
 71   explicit G4VMscModel(const G4String& nam);
 72 
 73   ~G4VMscModel() override;
 74 
 75   virtual G4double ComputeTruePathLengthLimit(const G4Track& track,  
 76                 G4double& stepLimit) = 0;
 77 
 78   virtual G4double ComputeGeomPathLength(G4double truePathLength) = 0;
 79 
 80   virtual G4double ComputeTrueStepLength(G4double geomPathLength) = 0;
 81 
 82   virtual G4ThreeVector& SampleScattering(const G4ThreeVector&,
 83             G4double safety) = 0;
 84 
 85   void InitialiseParameters(const G4ParticleDefinition*);
 86 
 87   void DumpParameters(std::ostream& out) const;
 88 
 89   // empty method
 90   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
 91        const G4MaterialCutsCouple*,
 92        const G4DynamicParticle*,
 93        G4double tmin, G4double tmax) override;
 94 
 95   //================================================================
 96   //  Set parameters of multiple scattering models
 97   //================================================================
 98  
 99   inline void SetStepLimitType(G4MscStepLimitType);
100 
101   inline void SetLateralDisplasmentFlag(G4bool val);
102 
103   inline void SetRangeFactor(G4double);
104 
105   inline void SetGeomFactor(G4double);
106 
107   inline void SetSkin(G4double);
108 
109   inline void SetLambdaLimit(G4double);
110 
111   inline void SetSafetyFactor(G4double);
112 
113   inline void SetSampleZ(G4bool);
114 
115   //================================================================
116   //  Get/Set access to Physics Tables
117   //================================================================
118 
119   inline G4VEnergyLossProcess* GetIonisation() const;
120 
121   inline void SetIonisation(G4VEnergyLossProcess*, 
122           const G4ParticleDefinition* part);
123 
124   //================================================================
125   //  Run time methods
126   //================================================================
127 
128 protected:
129 
130   // initialisation of the ParticleChange for the model
131   // initialisation of interface with geometry and ionisation 
132   G4ParticleChangeForMSC* 
133   GetParticleChangeForMSC(const G4ParticleDefinition* p = nullptr);
134 
135   // convert true length to geometry length
136   inline G4double ConvertTrueToGeom(G4double& tLength, G4double& gLength);
137 
138 public:
139 
140   // compute safety
141   inline G4double ComputeSafety(const G4ThreeVector& position, 
142         G4double limit= DBL_MAX);
143 
144   // compute linear distance to a geometry boundary
145   inline G4double ComputeGeomLimit(const G4Track&, G4double& presafety, 
146            G4double limit);
147 
148   inline G4double GetDEDX(const G4ParticleDefinition* part,
149                           G4double kineticEnergy,
150                           const G4MaterialCutsCouple* couple);
151 
152   inline G4double GetDEDX(const G4ParticleDefinition* part,
153                           G4double kineticEnergy,
154                           const G4MaterialCutsCouple* couple,
155                           G4double logKineticEnergy);
156 
157   inline G4double GetRange(const G4ParticleDefinition* part,
158                            G4double kineticEnergy,
159                            const G4MaterialCutsCouple* couple);
160 
161   inline G4double GetRange(const G4ParticleDefinition* part,
162                            G4double kineticEnergy,
163                            const G4MaterialCutsCouple* couple,
164                            G4double logKineticEnergy);
165 
166   inline G4double GetEnergy(const G4ParticleDefinition* part,
167           G4double range,
168           const G4MaterialCutsCouple* couple);
169 
170   // G4MaterialCutsCouple should be defined before call to this method
171   inline 
172   G4double GetTransportMeanFreePath(const G4ParticleDefinition* part,
173                                     G4double kinEnergy);
174 
175   inline
176   G4double GetTransportMeanFreePath(const G4ParticleDefinition* part,
177                                     G4double kinEnergy,
178                                     G4double logKinEnergy);
179 
180   //  hide assignment operator
181   G4VMscModel & operator=(const  G4VMscModel &right) = delete;
182   G4VMscModel(const  G4VMscModel&) = delete;
183 
184 private:
185 
186   G4SafetyHelper* safetyHelper = nullptr;
187   G4VEnergyLossProcess* ionisation = nullptr;
188   const G4ParticleDefinition* currentPart = nullptr;
189 
190   G4double dedx = 0.0;
191   G4double localtkin = 0.0;
192   G4double localrange = DBL_MAX;
193 
194 protected:
195 
196   G4double facrange = 0.04;
197   G4double facgeom = 2.5;
198   G4double facsafety = 0.6;
199   G4double skin = 1.0;
200   G4double dtrl = 0.05;
201   G4double lambdalimit;
202   G4double geomMin;
203   G4double geomMax;
204 
205   G4bool   samplez = false;
206   G4bool   latDisplasment = true;
207 
208   G4ThreeVector      fDisplacement;
209   G4MscStepLimitType steppingAlgorithm;
210 
211 };
212 
213 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215 
216 inline void G4VMscModel::SetLateralDisplasmentFlag(G4bool val)
217 {
218   if(!IsLocked()) { latDisplasment = val; }
219 }
220 
221 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
222 
223 inline void G4VMscModel::SetSkin(G4double val)
224 {
225   if(!IsLocked()) { skin = val; }
226 }
227 
228 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
229 
230 inline void G4VMscModel::SetRangeFactor(G4double val)
231 {
232   if(!IsLocked()) { facrange = val; }
233 }
234 
235 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
236 
237 inline void G4VMscModel::SetGeomFactor(G4double val)
238 {
239   if(!IsLocked()) { facgeom = val; }
240 }
241 
242 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
243 
244 inline void G4VMscModel::SetLambdaLimit(G4double val)
245 {
246   if(!IsLocked()) { lambdalimit = val; }
247 }
248 
249 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
250 
251 inline void G4VMscModel::SetSafetyFactor(G4double val)
252 {
253   if(!IsLocked()) { facsafety = val; }
254 }
255 
256 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
257 
258 inline void G4VMscModel::SetStepLimitType(G4MscStepLimitType val)
259 {
260   if(!IsLocked()) { steppingAlgorithm = val; }
261 }
262 
263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
264 
265 inline void G4VMscModel::SetSampleZ(G4bool val)
266 {
267   if(!IsLocked()) { samplez = val; }
268 }
269 
270 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
271 
272 inline G4double G4VMscModel::ComputeSafety(const G4ThreeVector& position, 
273              G4double limit)
274 {
275    return safetyHelper->ComputeSafety(position, limit);
276 }
277 
278 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
279 
280 inline G4double G4VMscModel::ConvertTrueToGeom(G4double& tlength, 
281                  G4double& glength)
282 {
283   glength = ComputeGeomPathLength(tlength);
284   // should return true length 
285   return tlength;
286 }
287 
288 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
289 
290 inline G4double G4VMscModel::ComputeGeomLimit(const G4Track& track, 
291                 G4double& presafety, 
292                 G4double limit)
293 {
294   return safetyHelper->CheckNextStep(
295           track.GetStep()->GetPreStepPoint()->GetPosition(),
296     track.GetMomentumDirection(),
297     limit, presafety);
298 }
299 
300 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
301 
302 inline G4double 
303 G4VMscModel::GetDEDX(const G4ParticleDefinition* part, G4double kinEnergy,
304                      const G4MaterialCutsCouple* couple)
305 {
306   G4double x;
307   if (nullptr != ionisation) {
308     x = ionisation->GetDEDX(kinEnergy, couple);
309   } else {
310     const G4double q = part->GetPDGCharge()*inveplus;
311     x = dedx*q*q;
312   }
313   return x;
314 }
315 
316 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
317 
318 inline G4double 
319 G4VMscModel::GetDEDX(const G4ParticleDefinition* part, G4double kinEnergy,
320                      const G4MaterialCutsCouple* couple, G4double logKinEnergy)
321 {
322   G4double x;
323   if (nullptr != ionisation) {
324     x = ionisation->GetDEDX(kinEnergy, couple, logKinEnergy);
325   } else {
326     const G4double q = part->GetPDGCharge()*inveplus;
327     x = dedx*q*q;
328   }
329   return x;
330 }
331 
332 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
333 
334 inline G4double 
335 G4VMscModel::GetRange(const G4ParticleDefinition* part, G4double kinEnergy,
336                       const G4MaterialCutsCouple* couple)
337 {
338   //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << "  " 
339   //  << ionisation << "  " << part->GetParticleName()
340   //  << G4endl;
341   localtkin = kinEnergy;
342   if (nullptr != ionisation) {
343     localrange = ionisation->GetRange(kinEnergy, couple); 
344   } else {
345     const G4double q = part->GetPDGCharge()*inveplus;
346     localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity()); 
347   }
348   //G4cout << "R(mm)= " << localrange << "  "  << ionisation << G4endl;
349   return localrange;
350 }
351 
352 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
353 
354 inline G4double 
355 G4VMscModel::GetRange(const G4ParticleDefinition* part,G4double kinEnergy, 
356                       const G4MaterialCutsCouple* couple, G4double logKinEnergy)
357 {
358   //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << "  " 
359   //  << ionisation << "  " << part->GetParticleName()
360   //  << G4endl;
361   localtkin = kinEnergy;
362   if (nullptr != ionisation) { 
363     localrange = ionisation->GetRange(kinEnergy, couple, logKinEnergy);
364   } else { 
365     const G4double q = part->GetPDGCharge()*inveplus;
366     localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
367   }
368   //G4cout << "R(mm)= " << localrange << "  "  << ionisation << G4endl;
369   return localrange;
370 }
371 
372 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
373 
374 inline G4double 
375 G4VMscModel::GetEnergy(const G4ParticleDefinition* part,
376            G4double range, const G4MaterialCutsCouple* couple)
377 {
378   G4double e;
379   //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << "  " << ionisation
380   //   << "  Rlocal(mm)= " << localrange << "  Elocal(MeV)= " << localtkin
381   //   << G4endl;
382   if(nullptr != ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
383   else { 
384     e = localtkin;
385     if(localrange > range) {
386       G4double q = part->GetPDGCharge()*inveplus;
387       e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity(); 
388     } 
389   }
390   return e;
391 }
392 
393 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
394 
395 inline G4VEnergyLossProcess* G4VMscModel::GetIonisation() const
396 {
397   return ionisation;
398 }
399 
400 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
401 
402 inline void G4VMscModel::SetIonisation(G4VEnergyLossProcess* p,
403                const G4ParticleDefinition* part)
404 {
405   ionisation = p;
406   currentPart = part;
407 }
408 
409 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
410 
411 inline G4double 
412 G4VMscModel::GetTransportMeanFreePath(const G4ParticleDefinition* part,
413                                       G4double ekin)
414 {
415   G4double x;
416   if (nullptr != xSectionTable) {
417     x = pFactor*(*xSectionTable)[basedCoupleIndex]->Value(ekin)/(ekin*ekin);
418   } else { 
419     x = pFactor*CrossSectionPerVolume(pBaseMaterial, part, ekin, 0.0, DBL_MAX); 
420   }
421   return (x > 0.0) ? 1.0/x : DBL_MAX;
422 }
423 
424 inline G4double 
425 G4VMscModel::GetTransportMeanFreePath(const G4ParticleDefinition* part,
426                                       G4double ekin, G4double logekin)
427 {
428   G4double x;
429   if (nullptr != xSectionTable) {
430     x = pFactor*(*xSectionTable)[basedCoupleIndex]->LogVectorValue(ekin, logekin)/(ekin*ekin);
431   } else { 
432     x = pFactor*CrossSectionPerVolume(pBaseMaterial, part, ekin, 0.0, DBL_MAX);
433   }
434   return (x > 0.0) ? 1.0/x : DBL_MAX;
435 }
436 
437 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
438 
439 #endif
440