Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/Geant4/G4VMscModel.hh was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // -------------------------------------------------------------------
0027 //
0028 // GEANT4 Class header file
0029 //
0030 //
0031 // File name:     G4VMscModel
0032 //
0033 // Author:        Vladimir Ivanchenko
0034 //
0035 // Creation date: 07.03.2008
0036 //
0037 // Modifications:
0038 // 07.04.2009 V.Ivanchenko moved msc methods from G4VEmModel to G4VMscModel 
0039 // 26.03.2012 V.Ivanchenko added transport x-section pointer and Get?Set methods
0040 //
0041 // Class Description:
0042 //
0043 // General interface to msc models
0044 
0045 // -------------------------------------------------------------------
0046 //
0047 #ifndef G4VMscModel_h
0048 #define G4VMscModel_h 1
0049 
0050 #include <CLHEP/Units/SystemOfUnits.h>
0051 
0052 #include "G4VEmModel.hh"
0053 #include "G4MscStepLimitType.hh"
0054 #include "globals.hh"
0055 #include "G4ThreeVector.hh"
0056 #include "G4Track.hh"
0057 #include "G4SafetyHelper.hh"
0058 #include "G4PhysicsTable.hh"
0059 #include "G4ThreeVector.hh"
0060 #include <vector>
0061 
0062 class G4ParticleChangeForMSC;
0063 class G4ParticleDefinition;
0064 class G4VEnergyLossProcess;
0065 
0066 class G4VMscModel : public G4VEmModel
0067 {
0068 
0069 public:
0070 
0071   explicit G4VMscModel(const G4String& nam);
0072 
0073   ~G4VMscModel() override;
0074 
0075   virtual G4double ComputeTruePathLengthLimit(const G4Track& track,  
0076                           G4double& stepLimit) = 0;
0077 
0078   virtual G4double ComputeGeomPathLength(G4double truePathLength) = 0;
0079 
0080   virtual G4double ComputeTrueStepLength(G4double geomPathLength) = 0;
0081 
0082   virtual G4ThreeVector& SampleScattering(const G4ThreeVector&,
0083                       G4double safety) = 0;
0084 
0085   void InitialiseParameters(const G4ParticleDefinition*);
0086 
0087   void DumpParameters(std::ostream& out) const;
0088 
0089   // empty method
0090   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0091              const G4MaterialCutsCouple*,
0092              const G4DynamicParticle*,
0093              G4double tmin, G4double tmax) override;
0094 
0095   //================================================================
0096   //  Set parameters of multiple scattering models
0097   //================================================================
0098  
0099   inline void SetStepLimitType(G4MscStepLimitType);
0100 
0101   inline void SetLateralDisplasmentFlag(G4bool val);
0102 
0103   inline void SetRangeFactor(G4double);
0104 
0105   inline void SetGeomFactor(G4double);
0106 
0107   inline void SetSkin(G4double);
0108 
0109   inline void SetLambdaLimit(G4double);
0110 
0111   inline void SetSafetyFactor(G4double);
0112 
0113   inline void SetSampleZ(G4bool);
0114 
0115   //================================================================
0116   //  Get/Set access to Physics Tables
0117   //================================================================
0118 
0119   inline G4VEnergyLossProcess* GetIonisation() const;
0120 
0121   inline void SetIonisation(G4VEnergyLossProcess*, 
0122                 const G4ParticleDefinition* part);
0123 
0124   //================================================================
0125   //  Run time methods
0126   //================================================================
0127 
0128 protected:
0129 
0130   // initialisation of the ParticleChange for the model
0131   // initialisation of interface with geometry and ionisation 
0132   G4ParticleChangeForMSC* 
0133   GetParticleChangeForMSC(const G4ParticleDefinition* p = nullptr);
0134 
0135   // convert true length to geometry length
0136   inline G4double ConvertTrueToGeom(G4double& tLength, G4double& gLength);
0137 
0138   // should be set before initialisation
0139   inline void SetUseSplineForMSC(G4bool val);
0140 
0141 public:
0142 
0143   // compute safety
0144   inline G4double ComputeSafety(const G4ThreeVector& position, 
0145                 G4double limit= DBL_MAX);
0146 
0147   // compute linear distance to a geometry boundary
0148   inline G4double ComputeGeomLimit(const G4Track&, G4double& presafety, 
0149                    G4double limit);
0150 
0151   G4double GetDEDX(const G4ParticleDefinition* part,
0152                           G4double kineticEnergy,
0153                           const G4MaterialCutsCouple* couple);
0154 
0155   G4double GetDEDX(const G4ParticleDefinition* part,
0156                           G4double kineticEnergy,
0157                           const G4MaterialCutsCouple* couple,
0158                           G4double logKineticEnergy);
0159 
0160   G4double GetRange(const G4ParticleDefinition* part,
0161                            G4double kineticEnergy,
0162                            const G4MaterialCutsCouple* couple);
0163 
0164   G4double GetRange(const G4ParticleDefinition* part,
0165                            G4double kineticEnergy,
0166                            const G4MaterialCutsCouple* couple,
0167                            G4double logKineticEnergy);
0168 
0169   G4double GetEnergy(const G4ParticleDefinition* part,
0170                 G4double range,
0171                 const G4MaterialCutsCouple* couple);
0172 
0173   // G4MaterialCutsCouple should be defined before call to this method
0174   inline 
0175   G4double GetTransportMeanFreePath(const G4ParticleDefinition* part,
0176                                     G4double kinEnergy);
0177 
0178   inline
0179   G4double GetTransportMeanFreePath(const G4ParticleDefinition* part,
0180                                     G4double kinEnergy,
0181                                     G4double logKinEnergy);
0182 
0183   //  hide assignment operator
0184   G4VMscModel & operator=(const  G4VMscModel &right) = delete;
0185   G4VMscModel(const  G4VMscModel&) = delete;
0186 
0187 private:
0188 
0189   G4SafetyHelper* safetyHelper = nullptr;
0190   G4VEnergyLossProcess* ionisation = nullptr;
0191   const G4ParticleDefinition* currentPart = nullptr;
0192 
0193   G4double dedx = 0.0;
0194   G4double localtkin = 0.0;
0195   G4double localrange = DBL_MAX;
0196 
0197 protected:
0198 
0199   G4double facrange = 0.04;
0200   G4double facgeom = 2.5;
0201   G4double facsafety = 0.6;
0202   G4double skin = 1.0;
0203   G4double dtrl = 0.05;
0204   G4double lambdalimit;
0205   G4double geomMin;
0206   G4double geomMax;
0207 
0208   G4ThreeVector fDisplacement;
0209   G4MscStepLimitType steppingAlgorithm;
0210 
0211   G4bool samplez = false;
0212   G4bool latDisplasment = true;
0213 
0214 private:
0215 
0216   G4bool useSpline = true;
0217 };
0218 
0219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0220 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0221 
0222 inline void G4VMscModel::SetLateralDisplasmentFlag(G4bool val)
0223 {
0224   if(!IsLocked()) { latDisplasment = val; }
0225 }
0226 
0227 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0228 
0229 inline void G4VMscModel::SetSkin(G4double val)
0230 {
0231   if(!IsLocked()) { skin = val; }
0232 }
0233 
0234 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0235 
0236 inline void G4VMscModel::SetRangeFactor(G4double val)
0237 {
0238   if(!IsLocked()) { facrange = val; }
0239 }
0240 
0241 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0242 
0243 inline void G4VMscModel::SetGeomFactor(G4double val)
0244 {
0245   if(!IsLocked()) { facgeom = val; }
0246 }
0247 
0248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0249 
0250 inline void G4VMscModel::SetLambdaLimit(G4double val)
0251 {
0252   if(!IsLocked()) { lambdalimit = val; }
0253 }
0254 
0255 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0256 
0257 inline void G4VMscModel::SetSafetyFactor(G4double val)
0258 {
0259   if(!IsLocked()) { facsafety = val; }
0260 }
0261 
0262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0263 
0264 inline void G4VMscModel::SetStepLimitType(G4MscStepLimitType val)
0265 {
0266   if(!IsLocked()) { steppingAlgorithm = val; }
0267 }
0268 
0269 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0270 
0271 inline void G4VMscModel::SetSampleZ(G4bool val)
0272 {
0273   if(!IsLocked()) { samplez = val; }
0274 }
0275 
0276 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0277 
0278 inline G4double G4VMscModel::ComputeSafety(const G4ThreeVector& position, 
0279                        G4double limit)
0280 {
0281    return safetyHelper->ComputeSafety(position, limit);
0282 }
0283 
0284 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0285 
0286 inline G4double G4VMscModel::ConvertTrueToGeom(G4double& tlength, 
0287                            G4double& glength)
0288 {
0289   glength = ComputeGeomPathLength(tlength);
0290   // should return true length 
0291   return tlength;
0292 }
0293 
0294 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0295 
0296 inline G4double G4VMscModel::ComputeGeomLimit(const G4Track& track, 
0297                           G4double& presafety, 
0298                           G4double limit)
0299 {
0300   return safetyHelper->CheckNextStep(
0301           track.GetStep()->GetPreStepPoint()->GetPosition(),
0302       track.GetMomentumDirection(),
0303       limit, presafety);
0304 }
0305 
0306 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0307 
0308 inline G4VEnergyLossProcess* G4VMscModel::GetIonisation() const
0309 {
0310   return ionisation;
0311 }
0312 
0313 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0314 
0315 inline void G4VMscModel::SetIonisation(G4VEnergyLossProcess* p,
0316                        const G4ParticleDefinition* part)
0317 {
0318   ionisation = p;
0319   currentPart = part;
0320 }
0321 
0322 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0323 
0324 inline G4double 
0325 G4VMscModel::GetTransportMeanFreePath(const G4ParticleDefinition* part,
0326                                       G4double ekin)
0327 {
0328   G4double x;
0329   if (nullptr != xSectionTable) {
0330     x = pFactor*(*xSectionTable)[basedCoupleIndex]->Value(ekin)/(ekin*ekin);
0331   } else { 
0332     x = pFactor*CrossSectionPerVolume(pBaseMaterial, part, ekin, 0.0, DBL_MAX); 
0333   }
0334   return (x > 0.0) ? 1.0/x : DBL_MAX;
0335 }
0336 
0337 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0338 
0339 inline G4double 
0340 G4VMscModel::GetTransportMeanFreePath(const G4ParticleDefinition* part,
0341                                       G4double ekin, G4double logekin)
0342 {
0343   G4double x;
0344   if (nullptr != xSectionTable) {
0345     x = pFactor*(*xSectionTable)[basedCoupleIndex]->LogVectorValue(ekin, logekin)/(ekin*ekin);
0346   } else { 
0347     x = pFactor*CrossSectionPerVolume(pBaseMaterial, part, ekin, 0.0, DBL_MAX);
0348   }
0349   return (x > 0.0) ? 1.0/x : DBL_MAX;
0350 }
0351 
0352 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0353 
0354 inline void G4VMscModel::SetUseSplineForMSC(G4bool val)
0355 {
0356   useSpline = val;
0357 }
0358 
0359 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0360 
0361 #endif