Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:59:25

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 //
0029 // GEANT4 Class header file
0030 //
0031 //
0032 // File name:     G4VMultipleScattering
0033 //
0034 // Author:        Vladimir Ivanchenko on base of Laszlo Urban code
0035 //
0036 // Creation date: 12.03.2002
0037 //
0038 // Modifications:
0039 //
0040 // 16-07-03 Update GetRange interface (V.Ivanchenko)
0041 // 26-11-03 bugfix in AlongStepDoIt (L.Urban)
0042 // 25-05-04 add protection against case when range is less than steplimit (VI)
0043 // 27-08-04 Add InitialiseForRun method (V.Ivanchneko)
0044 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
0045 // 15-04-05 remove boundary flag (V.Ivanchenko)
0046 // 07-10-05 error in a protection in GetContinuousStepLimit corrected (L.Urban)
0047 // 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko)
0048 // 26-01-06 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
0049 // 17-02-06 Save table of transport cross sections not mfp (V.Ivanchenko)
0050 // 07-03-06 Move step limit calculation to model (V.Ivanchenko)
0051 // 13-05-06 Add method to access model by index (V.Ivanchenko)
0052 // 12-02-07 Add get/set skin (V.Ivanchenko)
0053 // 27-10-07 Virtual functions moved to source (V.Ivanchenko)
0054 // 15-07-08 Reorder class members for further multi-thread development (VI)
0055 // 07-04-09 Moved msc methods from G4VEmModel to G4VMscModel (VI) 
0056 //
0057 // Class Description:
0058 //
0059 // It is the generic process of multiple scattering it includes common
0060 // part of calculations for all charged particles
0061 
0062 // -------------------------------------------------------------------
0063 //
0064 
0065 #ifndef G4VMultipleScattering_h
0066 #define G4VMultipleScattering_h 1
0067 
0068 #include "G4VContinuousDiscreteProcess.hh"
0069 #include "globals.hh"
0070 #include "G4Material.hh"
0071 #include "G4ParticleChangeForMSC.hh"
0072 #include "G4Track.hh"
0073 #include "G4Step.hh"
0074 #include "G4EmModelManager.hh"
0075 #include "G4VMscModel.hh"
0076 #include "G4EmParameters.hh"
0077 #include "G4MscStepLimitType.hh"
0078 
0079 class G4ParticleDefinition;
0080 class G4VEnergyLossProcess;
0081 class G4LossTableManager;
0082 class G4SafetyHelper;
0083 
0084 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0085 
0086 class G4VMultipleScattering : public G4VContinuousDiscreteProcess
0087 {
0088 public:
0089 
0090   explicit G4VMultipleScattering(const G4String& name = "msc",
0091                                  G4ProcessType type = fElectromagnetic);
0092 
0093   ~G4VMultipleScattering() override;
0094 
0095   //------------------------------------------------------------------------
0096   // Virtual methods to be implemented for the concrete model
0097   //------------------------------------------------------------------------
0098 
0099   void ProcessDescription(std::ostream& outFile) const override;
0100 
0101   virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
0102 
0103   // Print out of generic class parameters
0104   void StreamInfo(std::ostream& outFile, const G4ParticleDefinition&,
0105                   G4bool rst = false) const;
0106 
0107 protected:
0108 
0109   virtual void StreamProcessInfo(std::ostream&) const {};
0110 
0111 public:
0112 
0113   //------------------------------------------------------------------------
0114   // Generic methods common to all ContinuousDiscrete processes
0115   //------------------------------------------------------------------------
0116 
0117   // Initialise for build of tables
0118   void PreparePhysicsTable(const G4ParticleDefinition&) override;
0119 
0120   // Build physics table during initialisation
0121   void BuildPhysicsTable(const G4ParticleDefinition&) override;
0122 
0123   // Store PhysicsTable in a file.
0124   // Return false in case of failure at I/O
0125   G4bool StorePhysicsTable(const G4ParticleDefinition*,
0126                            const G4String& directory,
0127                            G4bool ascii = false) override;
0128 
0129   // Retrieve Physics from a file.
0130   // (return true if the Physics Table can be build by using file)
0131   // (return false if the process has no functionality or in case of failure)
0132   // File name should is constructed as processName+particleName and the
0133   // should be placed under the directory specified by the argument.
0134   G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
0135                               const G4String& directory,
0136                               G4bool ascii) override;
0137 
0138   // This is called in the beginning of tracking for a new track
0139   void StartTracking(G4Track*) override;
0140 
0141   // The function overloads the corresponding function of the base
0142   // class.It limits the step near to boundaries only
0143   // and invokes the method GetMscContinuousStepLimit at every step.
0144   G4double AlongStepGetPhysicalInteractionLength(
0145                                         const G4Track&,
0146                                         G4double  previousStepSize,
0147                                         G4double  currentMinimalStep,
0148                                         G4double& currentSafety,
0149                                         G4GPILSelection* selection) override;
0150 
0151   // The function overloads the corresponding function of the base
0152   // class.
0153   G4double PostStepGetPhysicalInteractionLength(
0154                                       const G4Track&,
0155                                       G4double  previousStepSize,
0156                                       G4ForceCondition* condition) override;
0157 
0158   // Along step actions
0159   G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override;
0160 
0161   // This method does not used for tracking, it is intended only for tests
0162   G4double ContinuousStepLimit(const G4Track& track,
0163                                G4double previousStepSize,
0164                                G4double currentMinimalStep,
0165                                G4double& currentSafety);
0166 
0167   // hide assignment operator
0168   G4VMultipleScattering(G4VMultipleScattering &) = delete;
0169   G4VMultipleScattering & operator=(const G4VMultipleScattering &right) = delete;
0170 
0171   //------------------------------------------------------------------------
0172   // Specific methods to set, access, modify models
0173   //------------------------------------------------------------------------
0174 
0175   // Select model in run time
0176   inline G4VEmModel* SelectModel(G4double kinEnergy, size_t idx);
0177 
0178 public:
0179 
0180   // Add model for region, smaller value of order defines which
0181   // model will be selected for a given energy interval  
0182   void AddEmModel(G4int order, G4VMscModel*, const G4Region* region = nullptr);
0183 
0184   // Assign a model to a process local list, to enable the list in run time 
0185   // the derived process should execute AddEmModel(..) for all such models
0186   void SetEmModel(G4VMscModel*, G4int idx = 0);
0187   
0188   // return a model from the local list
0189   inline G4VMscModel* EmModel(size_t index = 0) const;
0190 
0191   // Access to run time models 
0192   inline G4int NumberOfModels() const;
0193 
0194   inline G4VMscModel* GetModelByIndex(G4int idx, G4bool ver = false) const;
0195 
0196   //------------------------------------------------------------------------
0197   // Get/Set parameters for simulation of multiple scattering
0198   //------------------------------------------------------------------------
0199 
0200   inline G4bool LateralDisplasmentFlag() const;
0201   
0202   inline G4double Skin() const;
0203   
0204   inline G4double RangeFactor() const;
0205   
0206   inline G4double GeomFactor() const;
0207  
0208   inline G4double PolarAngleLimit() const;
0209 
0210   inline G4bool UseBaseMaterial() const;
0211  
0212   inline G4MscStepLimitType StepLimitType() const;
0213   inline void SetStepLimitType(G4MscStepLimitType val);
0214 
0215   inline G4double LowestKinEnergy() const;
0216   inline void SetLowestKinEnergy(G4double val);
0217 
0218   inline const G4ParticleDefinition* FirstParticle() const;
0219 
0220   //------------------------------------------------------------------------
0221   // Run time methods
0222   //------------------------------------------------------------------------
0223 
0224 protected:
0225 
0226   // This method is not used for tracking, it returns mean free path value
0227   G4double GetMeanFreePath(const G4Track& track,
0228                            G4double,
0229                            G4ForceCondition* condition) override;
0230 
0231   // This method is not used for tracking, it returns step limit
0232   G4double GetContinuousStepLimit(const G4Track& track,
0233                                   G4double previousStepSize,
0234                                   G4double currentMinimalStep,
0235                                   G4double& currentSafety) override;
0236 
0237 private:
0238 
0239   // ======== Parameters of the class fixed at construction =========
0240 
0241   G4EmModelManager*           modelManager;
0242   G4LossTableManager*         emManager;
0243   G4EmParameters*             theParameters;  
0244 
0245   // ======== Parameters of the class fixed at initialisation =======
0246 
0247   G4SafetyHelper*             safetyHelper = nullptr;
0248   const G4ParticleDefinition* firstParticle = nullptr;
0249   const G4ParticleDefinition* currParticle = nullptr;
0250  
0251   std::vector<G4VMscModel*>   mscModels;
0252 
0253   G4double                    facrange = 0.04;
0254   G4double                    lowestKinEnergy;
0255 
0256   // ======== Cached values - may be state dependent ================
0257 
0258 protected:
0259 
0260   G4ParticleChangeForMSC      fParticleChange;
0261 
0262 private:
0263 
0264   G4ThreeVector               fNewPosition;
0265   G4ThreeVector               fNewDirection;
0266 
0267   G4VMscModel*                currentModel = nullptr;
0268   G4VEnergyLossProcess*       fIonisation = nullptr;
0269 
0270   G4double                    geomMin;
0271   G4double                    minDisplacement2;
0272   G4double                    physStepLimit = 0.0;
0273   G4double                    tPathLength = 0.0;
0274   G4double                    gPathLength = 0.0;
0275 
0276   G4MscStepLimitType          stepLimit = fUseSafety;
0277   G4int                       numberOfModels = 0;
0278 
0279   G4bool                      latDisplacement = true;
0280   G4bool                      isIon = false;
0281   G4bool                      fPositionChanged = false;
0282   G4bool                      isActive = false;
0283   G4bool                      baseMat = false;
0284 };
0285 
0286 // ======== Run time inline methods ================
0287 
0288 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0289 
0290 inline G4VEmModel* 
0291 G4VMultipleScattering::SelectModel(G4double kinEnergy, size_t coupleIndex)
0292 {
0293   return modelManager->SelectModel(kinEnergy, coupleIndex);
0294 }
0295 
0296 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0297 
0298 inline  G4bool G4VMultipleScattering::LateralDisplasmentFlag() const
0299 {
0300   return latDisplacement;
0301 }
0302 
0303 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0304 
0305 inline  G4double G4VMultipleScattering::Skin() const
0306 {
0307   return theParameters->MscSkin();
0308 }
0309 
0310 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0311 
0312 inline  G4double G4VMultipleScattering::RangeFactor() const
0313 {
0314   return facrange;
0315 }
0316 
0317 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0318 
0319 inline  G4double G4VMultipleScattering::GeomFactor() const
0320 {
0321   return theParameters->MscGeomFactor();
0322 }
0323 
0324 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0325 
0326 inline  G4double G4VMultipleScattering::PolarAngleLimit() const
0327 {
0328   return theParameters->MscThetaLimit();
0329 }
0330 
0331 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0332 
0333 inline G4MscStepLimitType G4VMultipleScattering::StepLimitType() const
0334 {
0335   return stepLimit;
0336 }
0337 
0338 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0339 
0340 inline void G4VMultipleScattering::SetStepLimitType(G4MscStepLimitType val) 
0341 {
0342   theParameters->SetMscStepLimitType(val);
0343 }
0344 
0345 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0346 
0347 inline G4double G4VMultipleScattering::LowestKinEnergy() const
0348 {
0349   return lowestKinEnergy;
0350 }
0351 
0352 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0353 
0354 inline void G4VMultipleScattering::SetLowestKinEnergy(G4double val)
0355 {
0356   lowestKinEnergy = val;
0357 }
0358 
0359 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0360 
0361 inline const G4ParticleDefinition* G4VMultipleScattering::FirstParticle() const
0362 {
0363   return firstParticle;
0364 }
0365 
0366 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0367 
0368 inline G4VMscModel* G4VMultipleScattering::EmModel(size_t index) const
0369 {
0370   return (index < mscModels.size()) ? mscModels[index] : nullptr;
0371 }
0372 
0373 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0374 
0375 inline G4int G4VMultipleScattering::NumberOfModels() const
0376 {
0377   return numberOfModels;
0378 }
0379 
0380 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0381 
0382 inline G4VMscModel* 
0383 G4VMultipleScattering::GetModelByIndex(G4int idx, G4bool ver) const
0384 {
0385   // static cast is possible inside this class
0386   return static_cast<G4VMscModel*>(modelManager->GetModel(idx, ver));
0387 }
0388 
0389 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0390 
0391 inline G4bool G4VMultipleScattering::UseBaseMaterial() const
0392 {
0393   return baseMat;
0394 }
0395 
0396 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0397 
0398 #endif