[966] | 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 | // |
---|
[1055] | 26 | // $Id: G4VMscModel.hh,v 1.9 2009/04/07 18:39:47 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
[966] | 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class header file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4VMscModel |
---|
| 35 | // |
---|
| 36 | // Author: Vladimir Ivanchenko |
---|
| 37 | // |
---|
| 38 | // Creation date: 07.03.2008 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
[1055] | 41 | // 07.04.2009 V.Ivanchenko moved msc methods from G4VEmModel to G4VMscModel |
---|
[966] | 42 | // |
---|
| 43 | // Class Description: |
---|
| 44 | // |
---|
| 45 | // General interface to msc models |
---|
| 46 | |
---|
| 47 | // ------------------------------------------------------------------- |
---|
| 48 | // |
---|
| 49 | |
---|
| 50 | #ifndef G4VMscModel_h |
---|
| 51 | #define G4VMscModel_h 1 |
---|
| 52 | |
---|
| 53 | #include "G4VEmModel.hh" |
---|
| 54 | #include "G4MscStepLimitType.hh" |
---|
| 55 | #include "globals.hh" |
---|
[1055] | 56 | #include "G4ThreeVector.hh" |
---|
| 57 | #include "G4Track.hh" |
---|
| 58 | #include "G4SafetyHelper.hh" |
---|
[966] | 59 | |
---|
[1055] | 60 | class G4ParticleChangeForMSC; |
---|
| 61 | |
---|
[966] | 62 | class G4VMscModel : public G4VEmModel |
---|
| 63 | { |
---|
| 64 | |
---|
| 65 | public: |
---|
| 66 | |
---|
| 67 | G4VMscModel(const G4String& nam); |
---|
| 68 | |
---|
| 69 | virtual ~G4VMscModel(); |
---|
| 70 | |
---|
[1055] | 71 | virtual G4double ComputeTruePathLengthLimit(const G4Track& track, |
---|
| 72 | G4PhysicsTable* theLambdaTable, |
---|
| 73 | G4double currentMinimalStep); |
---|
| 74 | |
---|
| 75 | virtual G4double ComputeGeomPathLength(G4double truePathLength); |
---|
| 76 | |
---|
| 77 | virtual G4double ComputeTrueStepLength(G4double geomPathLength); |
---|
| 78 | |
---|
| 79 | virtual void SampleScattering(const G4DynamicParticle*, |
---|
| 80 | G4double safety); |
---|
| 81 | |
---|
| 82 | // empty |
---|
| 83 | virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*, |
---|
| 84 | const G4MaterialCutsCouple*, |
---|
| 85 | const G4DynamicParticle*, |
---|
| 86 | G4double tmin, |
---|
| 87 | G4double tmax); |
---|
| 88 | |
---|
| 89 | //================================================================ |
---|
| 90 | // Set parameters of multiple scattering models |
---|
| 91 | //================================================================ |
---|
| 92 | |
---|
[966] | 93 | inline void SetStepLimitType(G4MscStepLimitType); |
---|
| 94 | |
---|
| 95 | inline void SetLateralDisplasmentFlag(G4bool val); |
---|
| 96 | |
---|
| 97 | inline void SetRangeFactor(G4double); |
---|
| 98 | |
---|
| 99 | inline void SetGeomFactor(G4double); |
---|
| 100 | |
---|
| 101 | inline void SetSkin(G4double); |
---|
| 102 | |
---|
[1055] | 103 | inline void SetSampleZ(G4bool); |
---|
| 104 | |
---|
| 105 | protected: |
---|
| 106 | |
---|
| 107 | // initialisation of the ParticleChange for the model |
---|
| 108 | G4ParticleChangeForMSC* GetParticleChangeForMSC(); |
---|
| 109 | |
---|
| 110 | // initialisation of interface with geometry |
---|
| 111 | void InitialiseSafetyHelper(); |
---|
| 112 | |
---|
| 113 | // shift point of the track PostStep |
---|
| 114 | void ComputeDisplacement(G4ParticleChangeForMSC*, |
---|
| 115 | const G4ThreeVector& displDir, |
---|
| 116 | G4double displacement, |
---|
| 117 | G4double postsafety); |
---|
| 118 | |
---|
| 119 | // compute safety |
---|
| 120 | inline G4double ComputeSafety(const G4ThreeVector& position, G4double limit); |
---|
| 121 | |
---|
| 122 | // compute linear distance to a geometry boundary |
---|
| 123 | inline G4double ComputeGeomLimit(const G4Track& position, G4double& presafety, |
---|
| 124 | G4double limit); |
---|
| 125 | |
---|
[966] | 126 | private: |
---|
| 127 | |
---|
| 128 | // hide assignment operator |
---|
| 129 | G4VMscModel & operator=(const G4VMscModel &right); |
---|
| 130 | G4VMscModel(const G4VMscModel&); |
---|
| 131 | |
---|
[1055] | 132 | G4SafetyHelper* safetyHelper; |
---|
| 133 | |
---|
[966] | 134 | protected: |
---|
| 135 | |
---|
| 136 | G4double facrange; |
---|
| 137 | G4double facgeom; |
---|
| 138 | G4double facsafety; |
---|
| 139 | G4double skin; |
---|
| 140 | G4double dtrl; |
---|
| 141 | G4double lambdalimit; |
---|
[1055] | 142 | G4double geommax; |
---|
[966] | 143 | |
---|
| 144 | G4MscStepLimitType steppingAlgorithm; |
---|
| 145 | |
---|
| 146 | G4bool samplez; |
---|
| 147 | G4bool latDisplasment; |
---|
| 148 | |
---|
| 149 | }; |
---|
| 150 | |
---|
| 151 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 152 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 153 | |
---|
| 154 | inline void G4VMscModel::SetLateralDisplasmentFlag(G4bool val) |
---|
| 155 | { |
---|
| 156 | latDisplasment = val; |
---|
| 157 | } |
---|
| 158 | |
---|
| 159 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 160 | |
---|
| 161 | inline void G4VMscModel::SetSkin(G4double val) |
---|
| 162 | { |
---|
| 163 | skin = val; |
---|
| 164 | } |
---|
| 165 | |
---|
| 166 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 167 | |
---|
| 168 | inline void G4VMscModel::SetRangeFactor(G4double val) |
---|
| 169 | { |
---|
| 170 | facrange = val; |
---|
| 171 | } |
---|
| 172 | |
---|
| 173 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 174 | |
---|
| 175 | inline void G4VMscModel::SetGeomFactor(G4double val) |
---|
| 176 | { |
---|
| 177 | facgeom = val; |
---|
| 178 | } |
---|
| 179 | |
---|
| 180 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 181 | |
---|
| 182 | inline void G4VMscModel::SetStepLimitType(G4MscStepLimitType val) |
---|
| 183 | { |
---|
| 184 | steppingAlgorithm = val; |
---|
| 185 | } |
---|
| 186 | |
---|
| 187 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 188 | |
---|
[1055] | 189 | inline void G4VMscModel::SetSampleZ(G4bool val) |
---|
| 190 | { |
---|
| 191 | samplez = val; |
---|
| 192 | } |
---|
| 193 | |
---|
| 194 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 195 | |
---|
| 196 | inline G4double G4VMscModel::ComputeSafety(const G4ThreeVector& position, |
---|
| 197 | G4double) |
---|
| 198 | { |
---|
| 199 | return safetyHelper->ComputeSafety(position); |
---|
| 200 | } |
---|
| 201 | |
---|
| 202 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 203 | |
---|
| 204 | inline G4double G4VMscModel::ComputeGeomLimit(const G4Track& track, |
---|
| 205 | G4double& presafety, |
---|
| 206 | G4double limit) |
---|
| 207 | { |
---|
| 208 | G4double res = geommax; |
---|
| 209 | if(track.GetVolume() != safetyHelper->GetWorldVolume()) { |
---|
| 210 | res = safetyHelper->CheckNextStep( |
---|
| 211 | track.GetStep()->GetPreStepPoint()->GetPosition(), |
---|
| 212 | track.GetMomentumDirection(), |
---|
| 213 | limit, presafety); |
---|
| 214 | } |
---|
| 215 | return res; |
---|
| 216 | } |
---|
| 217 | |
---|
| 218 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 219 | |
---|
[966] | 220 | #endif |
---|
| 221 | |
---|