Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/Geant4/G4WentzelVIModel.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 //
0029 //
0030 // GEANT4 Class header file
0031 //
0032 //
0033 // File name:     G4WentzelVIModel
0034 //
0035 // Author:        V.Ivanchenko 
0036 //
0037 // Creation date: 09.04.2008 from G4MuMscModel
0038 //
0039 // Modifications:
0040 // 27-05-2010 V.Ivanchenko added G4WentzelOKandVIxSection class to
0041 //              compute cross sections and sample scattering angle
0042 //
0043 // Class Description:
0044 //
0045 // Implementation of the model of multiple scattering based on
0046 // G.Wentzel, Z. Phys. 40 (1927) 590.
0047 // H.W.Lewis, Phys Rev 78 (1950) 526.
0048 // J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
0049 // L.Urban, CERN-OPEN-2006-077.
0050 
0051 // -------------------------------------------------------------------
0052 //
0053 
0054 #ifndef G4WentzelVIModel_h
0055 #define G4WentzelVIModel_h 1
0056 
0057 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0058 
0059 #include "G4VMscModel.hh"
0060 #include "G4MaterialCutsCouple.hh"
0061 #include "G4WentzelOKandVIxSection.hh"
0062 
0063 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0064 
0065 class G4WentzelVIModel : public G4VMscModel
0066 {
0067 
0068 public:
0069 
0070   explicit G4WentzelVIModel(G4bool comb=true, const G4String& nam = "WentzelVIUni");
0071 
0072   ~G4WentzelVIModel() override;
0073 
0074   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0075 
0076   void InitialiseLocal(const G4ParticleDefinition*, 
0077                G4VEmModel* masterModel) override;
0078 
0079   void StartTracking(G4Track*) override;
0080 
0081   G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
0082                       G4double KineticEnergy,
0083                       G4double AtomicNumber,
0084                       G4double AtomicWeight=0., 
0085                       G4double cut = DBL_MAX,
0086                       G4double emax= DBL_MAX) override;
0087 
0088   G4ThreeVector& SampleScattering(const G4ThreeVector&, 
0089                   G4double safety) override;
0090 
0091   G4double 
0092   ComputeTruePathLengthLimit(const G4Track& track,
0093                  G4double& currentMinimalStep) override;
0094 
0095   G4double ComputeGeomPathLength(G4double truePathLength) override;
0096 
0097   G4double ComputeTrueStepLength(G4double geomStepLength) override;
0098 
0099   // defines low energy limit on energy transfer to atomic electron
0100   void SetFixedCut(G4double);
0101 
0102   // low energy limit on energy transfer to atomic electron
0103   G4double GetFixedCut() const;
0104 
0105   // access to cross section class
0106   void SetWVICrossSection(G4WentzelOKandVIxSection*);
0107 
0108   G4WentzelOKandVIxSection* GetWVICrossSection();
0109 
0110   void SetUseSecondMoment(G4bool);
0111 
0112   G4bool UseSecondMoment() const;
0113 
0114   G4PhysicsTable* GetSecondMomentTable();
0115 
0116   G4double SecondMoment(const G4ParticleDefinition*,
0117             const G4MaterialCutsCouple*,
0118             G4double kineticEnergy);
0119 
0120   void SetSingleScatteringFactor(G4double);
0121 
0122   void DefineMaterial(const G4MaterialCutsCouple*);
0123 
0124   G4WentzelVIModel & operator=(const G4WentzelVIModel &right) = delete;
0125   G4WentzelVIModel(const G4WentzelVIModel&) = delete;
0126 
0127 protected:
0128 
0129   G4double ComputeTransportXSectionPerVolume(G4double cosTheta);
0130 
0131   inline void SetupParticle(const G4ParticleDefinition*);
0132 
0133 private:
0134 
0135   G4double ComputeSecondMoment(const G4ParticleDefinition*,
0136                    G4double kineticEnergy);
0137 
0138 protected:
0139 
0140   G4WentzelOKandVIxSection* wokvi;
0141   const G4MaterialCutsCouple* currentCouple = nullptr;
0142   const G4Material* currentMaterial = nullptr;
0143 
0144   const G4ParticleDefinition* particle = nullptr;
0145   G4ParticleChangeForMSC* fParticleChange = nullptr;
0146   const G4DataVector* currentCuts = nullptr;
0147   G4PhysicsTable* fSecondMoments = nullptr;
0148 
0149   G4double lowEnergyLimit;
0150   G4double tlimitminfix;
0151   G4double ssFactor = 1.05;
0152   G4double invssFactor = 1.0;
0153 
0154   // cache kinematics
0155   G4double preKinEnergy = 0.0;
0156   G4double tPathLength = 0.0;
0157   G4double zPathLength = 0.0;
0158   G4double lambdaeff = 0.0;
0159   G4double currentRange = 0.0; 
0160   G4double cosTetMaxNuc = 0.0;
0161 
0162   G4double fixedCut = -1.0;
0163 
0164   // cache kinematics
0165   G4double effKinEnergy = 0.0;
0166 
0167   // single scattering parameters
0168   G4double cosThetaMin = 1.0;
0169   G4double cosThetaMax = -1.0;
0170   G4double xtsec = 0.0;
0171 
0172   G4int  currentMaterialIndex = 0;
0173   size_t idx2 = 0;
0174 
0175   // data for single scattering mode
0176   G4int nelments = 0;
0177 
0178   // flags
0179   G4bool   singleScatteringMode;
0180   G4bool   isCombined;
0181   G4bool   useSecondMoment;
0182 
0183   std::vector<G4double> xsecn;
0184   std::vector<G4double> prob;
0185 };
0186 
0187 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0189 
0190 inline void G4WentzelVIModel::SetupParticle(const G4ParticleDefinition* p)
0191 {
0192   // Initialise mass and charge
0193   if(p != particle) {
0194     particle = p;
0195     wokvi->SetupParticle(p);
0196   }
0197 }
0198 
0199 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0200 
0201 inline void G4WentzelVIModel::SetFixedCut(G4double val)
0202 {
0203   fixedCut = val;
0204 }
0205 
0206 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0207 
0208 inline G4double G4WentzelVIModel::GetFixedCut() const
0209 {
0210   return fixedCut;
0211 }
0212 
0213 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0214 
0215 inline void G4WentzelVIModel::SetWVICrossSection(G4WentzelOKandVIxSection* ptr)
0216 {
0217   if(ptr != wokvi) {
0218     delete wokvi;
0219     wokvi = ptr;
0220   }
0221 }
0222 
0223 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0224 
0225 inline G4WentzelOKandVIxSection* G4WentzelVIModel::GetWVICrossSection()
0226 {
0227   return wokvi;
0228 }
0229 
0230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0231 
0232 inline void G4WentzelVIModel::SetUseSecondMoment(G4bool val)
0233 {
0234   useSecondMoment = val;
0235 }
0236 
0237 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0238 
0239 inline G4bool G4WentzelVIModel::UseSecondMoment() const
0240 {
0241   return useSecondMoment;
0242 }
0243 
0244 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0245 
0246 inline G4PhysicsTable* G4WentzelVIModel::GetSecondMomentTable()
0247 {
0248   return fSecondMoments;
0249 }
0250 
0251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0252 
0253 inline G4double 
0254 G4WentzelVIModel::SecondMoment(const G4ParticleDefinition* part,
0255                    const G4MaterialCutsCouple* couple,
0256                    G4double ekin)
0257 {
0258   G4double x = 0.0;
0259   if(useSecondMoment) { 
0260     DefineMaterial(couple);
0261     x = (fSecondMoments) ?  
0262       (*fSecondMoments)[(*theDensityIdx)[currentMaterialIndex]]->Value(ekin, idx2)
0263       *(*theDensityFactor)[currentMaterialIndex]/(ekin*ekin)
0264       : ComputeSecondMoment(part, ekin);
0265   }
0266   return x;
0267 }
0268 
0269 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0270 
0271 #endif
0272