| [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 | //
|
|---|
| [1007] | 26 | // $Id: G4WentzelVIModel.hh,v 1.7 2008/08/04 08:49:09 vnivanch Exp $
|
|---|
| 27 | // GEANT4 tag $Name: geant4-09-02 $
|
|---|
| [966] | 28 | //
|
|---|
| 29 | // -------------------------------------------------------------------
|
|---|
| 30 | //
|
|---|
| 31 | //
|
|---|
| 32 | // GEANT4 Class header file
|
|---|
| 33 | //
|
|---|
| 34 | //
|
|---|
| 35 | // File name: G4WentzelVIModel
|
|---|
| 36 | //
|
|---|
| 37 | // Author: V.Ivanchenko
|
|---|
| 38 | //
|
|---|
| 39 | // Creation date: 09.04.2008 from G4MuMscModel
|
|---|
| 40 | //
|
|---|
| 41 | // Modifications:
|
|---|
| 42 | //
|
|---|
| 43 | //
|
|---|
| 44 | // Class Description:
|
|---|
| 45 | //
|
|---|
| 46 | // Implementation of the model of multiple scattering based on
|
|---|
| 47 | // G.Wentzel, Z. Phys. 40 (1927) 590.
|
|---|
| 48 | // H.W.Lewis, Phys Rev 78 (1950) 526.
|
|---|
| 49 | // J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
|
|---|
| 50 | // L.Urban, CERN-OPEN-2006-077.
|
|---|
| 51 |
|
|---|
| 52 | // -------------------------------------------------------------------
|
|---|
| 53 | //
|
|---|
| 54 |
|
|---|
| 55 | #ifndef G4WentzelVIModel_h
|
|---|
| 56 | #define G4WentzelVIModel_h 1
|
|---|
| 57 |
|
|---|
| 58 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 59 |
|
|---|
| 60 | #include "G4VMscModel.hh"
|
|---|
| 61 | #include "G4PhysicsTable.hh"
|
|---|
| 62 | #include "G4MscStepLimitType.hh"
|
|---|
| 63 | #include "G4MaterialCutsCouple.hh"
|
|---|
| 64 | #include "G4NistManager.hh"
|
|---|
| 65 |
|
|---|
| 66 | class G4LossTableManager;
|
|---|
| 67 | class G4ParticleChangeForMSC;
|
|---|
| [1007] | 68 | class G4SafetyHelper;
|
|---|
| [966] | 69 | class G4ParticleDefinition;
|
|---|
| 70 |
|
|---|
| 71 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 72 |
|
|---|
| 73 | class G4WentzelVIModel : public G4VMscModel
|
|---|
| 74 | {
|
|---|
| 75 |
|
|---|
| 76 | public:
|
|---|
| 77 |
|
|---|
| 78 | G4WentzelVIModel(const G4String& nam = "WentzelVIUni");
|
|---|
| 79 |
|
|---|
| 80 | virtual ~G4WentzelVIModel();
|
|---|
| 81 |
|
|---|
| 82 | void Initialise(const G4ParticleDefinition*, const G4DataVector&);
|
|---|
| 83 |
|
|---|
| 84 | G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
|
|---|
| 85 | G4double KineticEnergy,
|
|---|
| 86 | G4double AtomicNumber,
|
|---|
| 87 | G4double AtomicWeight=0.,
|
|---|
| 88 | G4double cut = DBL_MAX,
|
|---|
| 89 | G4double emax= DBL_MAX);
|
|---|
| 90 |
|
|---|
| 91 | void SampleScattering(const G4DynamicParticle*, G4double safety);
|
|---|
| 92 |
|
|---|
| 93 | void SampleSecondaries(std::vector<G4DynamicParticle*>*,
|
|---|
| 94 | const G4MaterialCutsCouple*,
|
|---|
| 95 | const G4DynamicParticle*,
|
|---|
| 96 | G4double,
|
|---|
| 97 | G4double);
|
|---|
| 98 |
|
|---|
| 99 | G4double ComputeTruePathLengthLimit(const G4Track& track,
|
|---|
| 100 | G4PhysicsTable* theLambdaTable,
|
|---|
| 101 | G4double currentMinimalStep);
|
|---|
| 102 |
|
|---|
| 103 | G4double ComputeGeomPathLength(G4double truePathLength);
|
|---|
| 104 |
|
|---|
| 105 | G4double ComputeTrueStepLength(G4double geomStepLength);
|
|---|
| 106 |
|
|---|
| 107 | private:
|
|---|
| 108 |
|
|---|
| 109 | G4double ComputeTransportXSectionPerVolume();
|
|---|
| 110 |
|
|---|
| 111 | G4double ComputeXSectionPerVolume();
|
|---|
| 112 |
|
|---|
| 113 | void ComputeMaxElectronScattering(G4double cut);
|
|---|
| 114 |
|
|---|
| 115 | inline G4double GetLambda(G4double kinEnergy);
|
|---|
| 116 |
|
|---|
| 117 | inline void SetupParticle(const G4ParticleDefinition*);
|
|---|
| 118 |
|
|---|
| 119 | inline void SetupKinematic(G4double kinEnergy, G4double cut);
|
|---|
| 120 |
|
|---|
| 121 | inline void SetupTarget(G4double Z, G4double kinEnergy);
|
|---|
| 122 |
|
|---|
| 123 | inline void DefineMaterial(const G4MaterialCutsCouple*);
|
|---|
| 124 |
|
|---|
| 125 | // hide assignment operator
|
|---|
| 126 | G4WentzelVIModel & operator=(const G4WentzelVIModel &right);
|
|---|
| 127 | G4WentzelVIModel(const G4WentzelVIModel&);
|
|---|
| 128 |
|
|---|
| 129 | const G4ParticleDefinition* theProton;
|
|---|
| 130 | const G4ParticleDefinition* theElectron;
|
|---|
| 131 | const G4ParticleDefinition* thePositron;
|
|---|
| 132 |
|
|---|
| 133 | G4ParticleChangeForMSC* fParticleChange;
|
|---|
| 134 |
|
|---|
| [1007] | 135 | G4SafetyHelper* safetyHelper;
|
|---|
| [966] | 136 | G4PhysicsTable* theLambdaTable;
|
|---|
| 137 | G4PhysicsTable* theLambda2Table;
|
|---|
| 138 | G4LossTableManager* theManager;
|
|---|
| 139 | const G4DataVector* currentCuts;
|
|---|
| 140 |
|
|---|
| 141 | G4NistManager* fNistManager;
|
|---|
| 142 |
|
|---|
| 143 | G4double numlimit;
|
|---|
| 144 | G4double tlimitminfix;
|
|---|
| 145 | G4double invsqrt12;
|
|---|
| 146 |
|
|---|
| 147 | // cash
|
|---|
| 148 | G4double preKinEnergy;
|
|---|
| 149 | G4double ecut;
|
|---|
| 150 | G4double lambda0;
|
|---|
| 151 | G4double tPathLength;
|
|---|
| 152 | G4double zPathLength;
|
|---|
| 153 | G4double lambdaeff;
|
|---|
| 154 | G4double currentRange;
|
|---|
| 155 | G4double par1;
|
|---|
| 156 | G4double par2;
|
|---|
| 157 | G4double par3;
|
|---|
| 158 |
|
|---|
| 159 | G4double xtsec;
|
|---|
| 160 | std::vector<G4double> xsecn;
|
|---|
| 161 | std::vector<G4double> prob;
|
|---|
| 162 | G4int nelments;
|
|---|
| 163 |
|
|---|
| 164 | G4int nbins;
|
|---|
| 165 | G4int nwarnings;
|
|---|
| 166 | G4int nwarnlimit;
|
|---|
| 167 |
|
|---|
| 168 | G4int currentMaterialIndex;
|
|---|
| 169 |
|
|---|
| 170 | const G4MaterialCutsCouple* currentCouple;
|
|---|
| 171 | const G4Material* currentMaterial;
|
|---|
| 172 |
|
|---|
| 173 | // single scattering parameters
|
|---|
| 174 | G4double coeff;
|
|---|
| 175 | G4double constn;
|
|---|
| 176 | G4double cosThetaMin;
|
|---|
| 177 | G4double cosThetaMax;
|
|---|
| 178 | G4double cosTetMaxNuc;
|
|---|
| 179 | G4double cosTetMaxNuc2;
|
|---|
| 180 | G4double cosTetMaxElec;
|
|---|
| 181 | G4double cosTetMaxElec2;
|
|---|
| 182 | G4double q2Limit;
|
|---|
| 183 | G4double alpha2;
|
|---|
| 184 | G4double a0;
|
|---|
| 185 |
|
|---|
| 186 | // projectile
|
|---|
| 187 | const G4ParticleDefinition* particle;
|
|---|
| 188 |
|
|---|
| 189 | G4double chargeSquare;
|
|---|
| 190 | G4double spin;
|
|---|
| 191 | G4double mass;
|
|---|
| 192 | G4double tkin;
|
|---|
| 193 | G4double mom2;
|
|---|
| 194 | G4double invbeta2;
|
|---|
| 195 | G4double etag;
|
|---|
| 196 | G4double lowEnergyLimit;
|
|---|
| 197 |
|
|---|
| 198 | // target
|
|---|
| 199 | G4double targetZ;
|
|---|
| 200 | G4double screenZ;
|
|---|
| 201 | G4double formfactA;
|
|---|
| 202 | G4double FF[100];
|
|---|
| 203 |
|
|---|
| 204 | // flags
|
|---|
| 205 | G4bool isInitialized;
|
|---|
| 206 | G4bool inside;
|
|---|
| 207 | };
|
|---|
| 208 |
|
|---|
| 209 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 210 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 211 |
|
|---|
| 212 | inline
|
|---|
| 213 | void G4WentzelVIModel::DefineMaterial(const G4MaterialCutsCouple* cup)
|
|---|
| 214 | {
|
|---|
| 215 | if(cup != currentCouple) {
|
|---|
| 216 | currentCouple = cup;
|
|---|
| 217 | currentMaterial = cup->GetMaterial();
|
|---|
| 218 | currentMaterialIndex = currentCouple->GetIndex();
|
|---|
| 219 | }
|
|---|
| 220 | }
|
|---|
| 221 |
|
|---|
| 222 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 223 |
|
|---|
| 224 | inline
|
|---|
| 225 | G4double G4WentzelVIModel::GetLambda(G4double e)
|
|---|
| 226 | {
|
|---|
| 227 | G4double x;
|
|---|
| 228 | if(theLambdaTable) {
|
|---|
| 229 | G4bool b;
|
|---|
| 230 | x = ((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b);
|
|---|
| 231 | } else {
|
|---|
| 232 | x = CrossSection(currentCouple,particle,e,
|
|---|
| 233 | (*currentCuts)[currentMaterialIndex]);
|
|---|
| 234 | }
|
|---|
| 235 | if(x > DBL_MIN) x = 1./x;
|
|---|
| 236 | else x = DBL_MAX;
|
|---|
| 237 | return x;
|
|---|
| 238 | }
|
|---|
| 239 |
|
|---|
| 240 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 241 |
|
|---|
| 242 | inline
|
|---|
| 243 | void G4WentzelVIModel::SetupParticle(const G4ParticleDefinition* p)
|
|---|
| 244 | {
|
|---|
| 245 | // Initialise mass and charge
|
|---|
| 246 | if(p != particle) {
|
|---|
| 247 | particle = p;
|
|---|
| 248 | mass = particle->GetPDGMass();
|
|---|
| 249 | spin = particle->GetPDGSpin();
|
|---|
| 250 | G4double q = particle->GetPDGCharge()/eplus;
|
|---|
| 251 | chargeSquare = q*q;
|
|---|
| 252 | tkin = 0.0;
|
|---|
| 253 | lowEnergyLimit = keV*mass/electron_mass_c2;
|
|---|
| 254 | }
|
|---|
| 255 | }
|
|---|
| 256 |
|
|---|
| 257 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 258 |
|
|---|
| 259 | inline void G4WentzelVIModel::SetupKinematic(G4double ekin, G4double cut)
|
|---|
| 260 | {
|
|---|
| 261 | if(ekin != tkin || ecut != cut) {
|
|---|
| 262 | tkin = ekin;
|
|---|
| 263 | mom2 = tkin*(tkin + 2.0*mass);
|
|---|
| 264 | invbeta2 = 1.0 + mass*mass/mom2;
|
|---|
| 265 | cosTetMaxNuc = cosThetaMax;
|
|---|
| 266 | if(ekin <= 10.*cut && mass < MeV) {
|
|---|
| 267 | cosTetMaxNuc = ekin*(cosThetaMax + 1.0)/(10.*cut) - 1.0;
|
|---|
| 268 | }
|
|---|
| 269 | ComputeMaxElectronScattering(cut);
|
|---|
| 270 | }
|
|---|
| 271 | }
|
|---|
| 272 |
|
|---|
| 273 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 274 |
|
|---|
| 275 | inline void G4WentzelVIModel::SetupTarget(G4double Z, G4double e)
|
|---|
| 276 | {
|
|---|
| 277 | if(Z != targetZ || e != etag) {
|
|---|
| 278 | etag = e;
|
|---|
| 279 | targetZ = Z;
|
|---|
| 280 | G4int iz= G4int(Z);
|
|---|
| 281 | if(iz > 99) iz = 99;
|
|---|
| 282 | G4double x = fNistManager->GetZ13(iz);
|
|---|
| 283 | screenZ = a0*x*x/mom2;
|
|---|
| 284 | if(iz > 1) screenZ *=(1.13 + 3.76*invbeta2*Z*Z*chargeSquare*alpha2);
|
|---|
| 285 | // screenZ = a0*x*x*(1.13 + 3.76*Z*Z*chargeSquare*alpha2)/mom2;
|
|---|
| 286 | // A.V. Butkevich et al., NIM A 488 (2002) 282
|
|---|
| 287 | formfactA = FF[iz];
|
|---|
| 288 | if(formfactA == 0.0) {
|
|---|
| 289 | x = fNistManager->GetA27(iz);
|
|---|
| 290 | formfactA = constn*x*x;
|
|---|
| 291 | FF[iz] = formfactA;
|
|---|
| 292 | }
|
|---|
| 293 | formfactA *= mom2;
|
|---|
| 294 | cosTetMaxNuc2 = cosTetMaxNuc;
|
|---|
| 295 | /*
|
|---|
| 296 | G4double ee = 10.*eV*Z;
|
|---|
| 297 | if(1 == iz) ee *= 2.0;
|
|---|
| 298 | G4double z = std::min(cosTetMaxElec, 1.0 - std::max(ecut,ee)*amu_c2
|
|---|
| 299 | *fNistManager->GetAtomicMassAmu(iz)/mom2);
|
|---|
| 300 | cosTetMaxElec2 = std::max(cosTetMaxNuc2, z);
|
|---|
| 301 | */
|
|---|
| 302 | cosTetMaxElec2 = cosTetMaxElec;
|
|---|
| 303 | }
|
|---|
| 304 | }
|
|---|
| 305 |
|
|---|
| 306 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 307 |
|
|---|
| 308 | #endif
|
|---|
| 309 |
|
|---|