Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4SpecialUrbanMscModel.h
Go to the documentation of this file.
1#ifndef TG4SpecialUrbanMscModel_h
2#define TG4SpecialUrbanMscModel_h 1
3
4//------------------------------------------------
5// The Geant4 Virtual Monte Carlo package
6// Copyright (C) 2007 - 2015 Ivana Hrivnacova
7// All rights reserved.
8//
9// For the licensing terms see geant4_vmc/LICENSE.
10// Contact: root-vmc@cern.ch
11//-------------------------------------------------
12
17
18//
19// ********************************************************************
20// * License and Disclaimer *
21// * *
22// * The Geant4 software is copyright of the Copyright Holders of *
23// * the Geant4 Collaboration. It is provided under the terms and *
24// * conditions of the Geant4 Software License, included in the file *
25// * LICENSE and available at http://cern.ch/geant4/license . These *
26// * include a list of copyright holders. *
27// * *
28// * Neither the authors of this software system, nor their employing *
29// * institutes,nor the agencies providing financial support for this *
30// * work make any representation or warranty, express or implied, *
31// * regarding this software system or assume any liability for its *
32// * use. Please see the license in the file LICENSE and URL above *
33// * for the full disclaimer and the limitation of liability. *
34// * *
35// * This code implementation is the result of the scientific and *
36// * technical work of the GEANT4 collaboration. *
37// * By using, copying, modifying or distributing the software (or *
38// * any work based on the software) you agree to acknowledge its *
39// * use in resulting scientific publications, and indicate your *
40// * acceptance of all terms of the Geant4 Software license. *
41// ********************************************************************
42//
43// $Id: $
44// GEANT4 tag $Name: $
45//
46// -------------------------------------------------------------------
47//
48//
49// GEANT4 Class header file
50//
51//
52// File name: TG4SpecialUrbanMscModel
53//
54// Author: Vladimir Ivanchenko adopt Laszlo Urban model for
55// ALICE requirements
56//
57// Creation date: 12.08.2015
58//
59
60// -------------------------------------------------------------------
61//
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
65#include <CLHEP/Units/SystemOfUnits.h>
66
67#include "G4Exp.hh"
68#include "G4Log.hh"
69#include "G4MscStepLimitType.hh"
70#include "G4VMscModel.hh"
71#include "G4Version.hh"
72
73class G4ParticleChangeForMSC;
74class G4SafetyHelper;
75class G4LossTableManager;
76
77namespace CLHEP
78{
79class HepRandomEngine;
80}
81
82static const G4double c_highland = 13.6 * CLHEP::MeV;
83
88
90{
91
92 public:
93 TG4SpecialUrbanMscModel(const G4String& nam = "ALICEUrbanMsc");
94
96
97 void Initialise(const G4ParticleDefinition*, const G4DataVector&);
98
99 void StartTracking(G4Track*);
100
102 G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight = 0.,
103 G4double cut = 0., G4double emax = DBL_MAX);
104
105 G4ThreeVector& SampleScattering(const G4ThreeVector&, G4double safety);
106
108 const G4Track& track, G4double& currentMinimalStep);
109
110 G4double ComputeGeomPathLength(G4double truePathLength);
111
112 G4double ComputeTrueStepLength(G4double geomStepLength);
113
114 inline G4double ComputeTheta0(
115 G4double truePathLength, G4double KineticEnergy);
116
117 private:
118 G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
119
120 inline void SetParticle(const G4ParticleDefinition*);
121
122 inline void UpdateCache();
123
124 inline G4double SimpleScattering(G4double xmeanth, G4double x2meanth);
125
126 inline G4double LatCorrelation();
127
128 // hide assignment operator
131
132#if G4VERSION_NUMBER >= 1020
133 CLHEP::HepRandomEngine* rndmEngineMod;
134#endif
135
137 G4ParticleChangeForMSC* fParticleChange;
138
139 const G4MaterialCutsCouple* couple;
140 G4LossTableManager* theManager;
141
142 G4double mass;
145
146 G4double taubig;
147 G4double tausmall;
148 G4double taulim;
149 G4double currentTau;
150 G4double tlimit;
151 G4double tlimitmin;
153 G4double tgeom;
154
155 G4double geombig;
156 G4double geommin;
157 G4double geomlimit;
158 G4double skindepth;
159 G4double smallstep;
160
161 G4double presafety;
162
163 G4double lambda0;
164 G4double lambdaeff;
165 G4double tPathLength;
166 G4double zPathLength;
167 G4double par1, par2, par3;
168
169 G4double stepmin;
170
172 G4double currentRange;
173 G4double rangeinit;
175
177 G4double third;
178
180
181 G4double y;
182 G4double Zold;
183 G4double Zeff, Z2, Z23, lnZ;
186
187 G4bool firstStep;
188 G4bool inside;
190
192
193 G4double rangecut;
194 G4double drr, finalr;
195};
196
197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
199
201{
202 if (p != particle) {
203 particle = p;
204 mass = p->GetPDGMass();
205 charge = p->GetPDGCharge() / CLHEP::eplus;
207 }
208}
209
210//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
211
213{
214 lnZ = G4Log(Zeff);
215 // correction in theta0 formula
216 G4double w = G4Exp(lnZ / 6.);
217 G4double facz = 0.990395 + w * (-0.168386 + w * 0.093286);
218 coeffth1 = facz * (1. - 8.7780e-2 / Zeff);
219 coeffth2 = facz * (4.0780e-2 + 1.7315e-4 * Zeff);
220
221 // tail parameters
222 G4double Z13 = w * w;
223 coeffc1 = 2.3785 - Z13 * (4.1981e-1 - Z13 * 6.3100e-2);
224 coeffc2 = 4.7526e-1 + Z13 * (1.7694 - Z13 * 3.3885e-1);
225 coeffc3 = 2.3683e-1 - Z13 * (1.8111 - Z13 * 3.2774e-1);
226 coeffc4 = 1.7888e-2 + Z13 * (1.9659e-2 - Z13 * 2.6664e-3);
227
228 Z2 = Zeff * Zeff;
229 Z23 = Z13 * Z13;
230
231 Zold = Zeff;
232}
233
234//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
235
237 G4double trueStepLength, G4double KineticEnergy)
238{
239 // for all particles take the width of the central part
240 // from a parametrization similar to the Highland formula
241 // ( Highland formula: Particle Physics Booklet, July 2002, eq. 26.10)
242 G4double invbetacp =
243 std::sqrt((currentKinEnergy + mass) * (KineticEnergy + mass) /
245 KineticEnergy * (KineticEnergy + 2. * mass)));
246 y = trueStepLength / currentRadLength;
247 G4double theta0 = c_highland * std::abs(charge) * std::sqrt(y) * invbetacp;
248 y = G4Log(y);
249 // correction factor from e- scattering data
250 theta0 *= (coeffth1 + coeffth2 * y);
251 return theta0;
252}
253
254//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
255
257 G4double xmeanth, G4double x2meanth)
258{
259 // 'large angle scattering'
260 // 2 model functions with correct xmean and x2mean
261 G4double a =
262 (2. * xmeanth + 9. * x2meanth - 3.) / (2. * xmeanth - 3. * x2meanth + 1.);
263 G4double prob = (a + 2.) * xmeanth / a;
264
265 // sampling
266 G4double cth = 1.;
267 if (rndmEngineMod->flat() < prob) {
268 cth = -1. + 2. * G4Exp(G4Log(rndmEngineMod->flat()) / (a + 1.));
269 }
270 else {
271 cth = -1. + 2. * rndmEngineMod->flat();
272 }
273 return cth;
274}
275
276//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
277
279{
280 static const G4double kappa = 2.5;
281 static const G4double kappami1 = 1.5;
282
283 G4double latcorr = 0.;
284 if ((currentTau >= tausmall) && !insideskin) {
285 if (currentTau < taulim)
286 latcorr = lambdaeff * kappa * currentTau * currentTau *
287 (1. - (kappa + 1.) * currentTau * third) * third;
288 else {
289 G4double etau = 0.;
290 if (currentTau < taubig) {
291 etau = G4Exp(-currentTau);
292 }
293 latcorr = -kappa * currentTau;
294 latcorr = G4Exp(latcorr) / kappami1;
295 latcorr += 1. - kappa * etau / kappami1;
296 latcorr *= 2. * lambdaeff * third;
297 }
298 }
299 return latcorr;
300}
301
302//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
303
304#endif
static const G4double c_highland
Laszlo Urban model adapted for ALICE EMCAL requirements.
void SetParticle(const G4ParticleDefinition *)
G4ParticleChangeForMSC * fParticleChange
TG4SpecialUrbanMscModel(const TG4SpecialUrbanMscModel &)
G4double ComputeGeomPathLength(G4double truePathLength)
const G4MaterialCutsCouple * couple
const G4ParticleDefinition * particle
G4double ComputeTrueStepLength(G4double geomStepLength)
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)
G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double SimpleScattering(G4double xmeanth, G4double x2meanth)
TG4SpecialUrbanMscModel & operator=(const TG4SpecialUrbanMscModel &right)
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX)
TG4SpecialUrbanMscModel(const G4String &nam="ALICEUrbanMsc")
G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy)