Garfield++ v2r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
HeedDeltaElectronCS.cpp
Go to the documentation of this file.
1#include <iomanip>
7
8// 2003, I. Smirnov
9
10namespace Heed {
11
12using CLHEP::twopi;
13using CLHEP::degree;
14using CLHEP::electron_mass_c2;
15using CLHEP::eV;
16using CLHEP::MeV;
17using CLHEP::cm;
18using CLHEP::cm3;
19using CLHEP::gram;
20using CLHEP::mole;
21using CLHEP::Avogadro;
22
24
26 ElElasticScat* fees,
28 PairProd* fpairprod, int fsruth,
29 double fmlambda, double fmthetac)
30 : hmd(fhmd),
31 ees(fees),
32 eesls(feesls),
33 pairprod(fpairprod),
34 mlambda(fmlambda),
35 sruth(fsruth),
36 mthetac(fmthetac) {
37 mfunname("HeedDeltaElectronCS::HeedDeltaElectronCS(...)");
38
39 const long qe = hmd->energy_mesh->get_q();
40 eLoss.resize(qe, 0.);
41 beta.resize(qe, 0.);
42 momentum.resize(qe, 0.);
43 std::vector<double> beta2(qe, 0.);
44 std::vector<double> momentum2(qe, 0.);
45
46 const double rho = hmd->matter->density();
47 const double zmean = hmd->matter->Z_mean();
48 const double amean = hmd->matter->A_mean();
49 const double ZA = zmean / amean;
50 const double I_eff = 15.8 * eV * zmean;
51 double smax = 0.0;
52 for (long ne = 0; ne < qe; ne++) {
53 const double ec = hmd->energy_mesh->get_ec(ne) * MeV;
54 const double gamma_1 = ec / electron_mass_c2;
55 beta[ne] = lorbeta(gamma_1);
56 beta2[ne] = beta[ne] * beta[ne];
57 const double en = electron_mass_c2 + ec;
58 momentum2[ne] =
59 (en * en - electron_mass_c2 * electron_mass_c2) / (MeV * MeV);
60 momentum[ne] = sqrt(momentum2[ne]);
61 const double dedx = e_cont_enloss(ZA, I_eff, rho, ec, DBL_MAX, -1);
62 if (smax < dedx) smax = dedx;
63 eLoss[ne] = dedx / (MeV / cm);
64 }
65 smax = smax / (MeV / cm);
66 // looking for minimal (maximal by value but negative) derivative
67 double deriv_min = 0.0;
68 long n_deriv_min = 0;
69 for (long ne = 1; ne < qe; ne++) {
70 double d =
71 (eLoss[ne] * beta[ne] - eLoss[ne - 1] * beta[ne - 1]) /
72 (hmd->energy_mesh->get_ec(ne) - hmd->energy_mesh->get_ec(ne - 1));
73 if (deriv_min > d) {
74 deriv_min = d;
75 n_deriv_min = ne - 1;
76 }
77 }
78 for (long ne = 0; ne < n_deriv_min - 1; ne++) {
79 eLoss[ne] = (eLoss[n_deriv_min - 1] * beta[n_deriv_min - 1] +
80 deriv_min * (hmd->energy_mesh->get_ec(ne) -
81 hmd->energy_mesh->get_ec(n_deriv_min - 1))) /
82 beta[ne];
83 }
84
85 lambda.resize(qe);
86 if (sruth == 2) {
88 angular_mesh_c[0] = 0.0;
89 const double amin = 0.3;
90 const double amax = 180.0;
91 const double rk = pow(amax / amin, (1.0 / (q_angular_mesh - 2.)));
92 double ar = amin;
93 angular_mesh_c[1] = ar;
94 for (long n = 2; n < q_angular_mesh; n++) {
95 angular_mesh_c[n] = angular_mesh_c[n - 1] * rk;
96 }
97 angular_mesh_c[q_angular_mesh - 1] = 180.0;
98 angular_points_ran.resize(qe);
99 low_angular_points_ran.resize(qe);
100 low_lambda.resize(qe);
101 const long qes = eesls->get_ees()->get_qe();
102#ifdef USE_MEAN_COEF
103 mean_coef_low_sigma.resize(qes);
104#else
105 coef_low_sigma.resize(qes);
106#endif
107 for (long ne = 0; ne < qes; ne++) {
108 long qat = hmd->matter->qatom();
109 double s = 0.0;
110 for (long nat = 0; nat < qat; nat++) {
111#ifdef USE_MEAN_COEF
112 s += eesls->get_mean_coef(hmd->matter->atom(nat)->Z(), ne) *
113 hmd->matter->weight_quan(nat);
114#else
115 s += eesls->get_coef(hmd->matter->atom(nat)->Z(), ne) *
116 hmd->matter->weight_quan(nat);
117#endif
118 }
119#ifdef USE_MEAN_COEF
120 mean_coef_low_sigma[ne] = s;
121#else
122 coef_low_sigma[ne] = s;
123#endif
124 }
125 }
126
127 for (long ne = 0; ne < qe; ne++) {
128 double rr;
129 double ek = hmd->energy_mesh->get_ec(ne) * 1000.0;
130 if (ek <= 10.0) {
131 rr = 1.0e-3 * amean / (gram / mole) / zmean * 3.872e-3 * pow(ek, 1.492);
132 } else {
133 rr = 1.0e-3 * 6.97e-3 * pow(ek, 1.6);
134 }
135 rr = rr / (rho / (gram / cm3));
136 rr = rr * 0.1;
137 // Iprintn(mcout, rr);
138 double cor = 1.0;
139 {
140 // b-k*(x-a)**2 = 0 => x= a +- sqrt(b/k)
141 // k = b / (x - a)**2
142 double a = 2.5;
143 double b = 4;
144 // k=1.0/4.0
145 double x = 0.0;
146 double k = b / ((x - a) * (x - a));
147 x = ek * 1000.0;
148 double r = b - k * (x - a) * (x - a);
149 if (r < 0.0)
150 r = 1;
151 else
152 r = r + 1;
153 cor = r;
154 }
155 if (sruth == 1) {
156 lambda[ne] = mlambda / (rho / (gram / cm3));
157 if (lambda[ne] < rr) lambda[ne] = rr;
158 lambda[ne] = lambda[ne] * cor;
159 // Calculate the minimum angle for restriction of field by atomic shell
160 double mT = 2.0 * asin(1.0 / (2.0 * momentum[ne] * zmean * 5.07e2));
161 // Throw out too slow interactions. They do not have any influence.
162 if (mT < mthetac) mT = mthetac;
163 // Calculate the cut angle due to mean free part
164 double A = hmd->Rutherford_const / cor / (momentum2[ne] * beta2[ne]) /
165 pow(5.07e10, 2);
166 double B = lambda[ne] * A;
167 B = sqrt(B / (B + 1.0));
168 // Threshold turn angle
169 double thetac = 2.0 * asin(B);
170
171 // If it too little, reset it. It will lead to increasing
172 // of lambda and decreasing of calculation time.
173 if (thetac < mT) {
174 thetac = mT;
175 B = mT; // B is double precision
176 double r = sin(0.5 * B);
177 lambda[ne] = 1 / A * 2.0 * r * r / (1 + cos(B));
178 // r=cos(TetacBdel(nen,nm))
179 // lamBdel=A*(1.0+r)/(1.0-r)
180 // lamBdel=1.0/lamBdel
181 // lamBdel=(p2*bet2*sin(TetacBdel/2.0)**2) / A
182 }
183 // const double CosThetac12 = cos(0.5 * thetac);
184 // const dobule SinThetac12 = sin(0.5 * thetac);
185 // Flag that scattering is spherical.
186 // const bool sisfera = thetac > 1.5 ? true : false;
187 } else if (sruth == 0) {
188 // Gauss formula
189 const double msig = sqrt(2.0) * 13.6 / (beta[ne] * momentum[ne]);
190 double x = mthetac / msig;
191 x = x * x;
192 x = x * hmd->radiation_length * cor;
193 lambda[ne] = mlambda / (rho / (gram / cm3));
194 if (lambda[ne] < rr) lambda[ne] = rr;
195 lambda[ne] = lambda[ne] * cor;
196 if (lambda[ne] < x) lambda[ne] = x;
197 } else if (sruth == 2) {
198 // Cross section for material per one av. atom, in cm^2/rad.
199 std::vector<double> smat(q_angular_mesh, 0.);
200 const double energy = hmd->energy_mesh->get_ec(ne);
201 const long qat = hmd->matter->qatom();
202 for (long nan = 0; nan < q_angular_mesh; nan++) {
203 const double angle = angular_mesh_c[nan] * degree;
204 double s = 0.0;
205 for (long nat = 0; nat < qat; nat++) {
206 const int z = hmd->matter->atom(nat)->Z();
207 const double w = hmd->matter->weight_quan(nat);
208 s += ees->get_CS(z, energy, angle) * w;
209 }
210 s = s * 1.0E-16;
211 // s = s * 1.0E-16 * C1_MEV_CM * C1_MEV_CM;
212 // Angstrem**2 -> cm**2
213 // cm**2 -> MeV**-2
214 smat[nan] = s * twopi * sin(angle); // sr -> dtheta
215 }
220 const double coef =
221 Avogadro * rho / (gram / cm3) / (amean / (gram / mole));
222 lambda[ne] =
223 1. / (angular_points_ran[ne].get_integ_active() * degree * coef);
224 low_lambda[ne] =
225 1. / (low_angular_points_ran[ne].get_integ_active() * degree * coef);
226 }
227 }
228}
229
230double HeedDeltaElectronCS::get_sigma(double energy, double nscat) const {
231 mfunname("double HeedDeltaElectronCS::get_sigma(...)");
232 check_econd11(nscat, < 0, mcerr);
233 // check_econd21(nscat , < 0 || , > eesls->get_qscat() , mcerr);
234 // ^ not compatible with Poisson
235 const long qe = ees->get_qe();
236 double energyKeV = energy * 1000.0;
237 energyKeV = std::max(energyKeV, ees->get_energy_mesh(0));
238 energyKeV = std::min(energyKeV, ees->get_energy_mesh(qe - 1));
239 long n1 = 0;
240 long n2 = qe - 1;
241 while (n2 - n1 > 1) {
242 const long n3 = n1 + (n2 - n1) / 2;
243 if (energyKeV < ees->get_energy_mesh(n3))
244 n2 = n3;
245 else
246 n1 = n3;
247 }
248 const double e1 = ees->get_energy_mesh(n1);
249 const double e2 = ees->get_energy_mesh(n2);
250#ifdef USE_MEAN_COEF
251 const double v1 = nscat * mean_coef_low_sigma[n1];
252 const double v2 = nscat * mean_coef_low_sigma[n2];
253#else
254 const double v1 = nscat * coef_low_sigma[n1];
255 const double v2 = nscat * coef_low_sigma[n2];
256#endif
257 return v1 + (v2 - v1) / (e2 - e1) * (energyKeV - e1);
258}
259
260void HeedDeltaElectronCS::print(std::ostream& file, int l) const {
261 if (l <= 0) return;
262 Ifile << "HeedDeltaElectronCS(l=" << l << "):";
263 long qe = hmd->energy_mesh->get_q();
264 // Iprintn(mcout, qe);
265 // mcout<<std::endl;
266 Iprintn(file, mlambda);
267 Iprintn(file, mthetac);
268 Iprintn(file, sruth);
269 Ifile << " get_ec, beta, momentum, eLoss, lambda, "
270 " low_lambda:" << std::endl;
271 indn.n += 2;
272 for (long ne = 0; ne < qe; ne++) {
273 Ifile << std::setw(3) << ne << ' ' << std::setw(12)
274 << hmd->energy_mesh->get_ec(ne) << ' ' << std::setw(12) << beta[ne]
275 << ' ' << std::setw(12) << momentum[ne] << ' ' << std::setw(12)
276 << eLoss[ne] << ' ' << std::setw(12) << lambda[ne] << ' '
277 << std::setw(12) << low_lambda[ne] << '\n';
278 }
279 indn.n -= 2;
280 Ifile << "na, angular_mesh_c:" << std::endl;
281 indn.n += 2;
282 for (long na = 0; na < q_angular_mesh; na++) {
283 Ifile << na << ' ' << std::setw(12) << angular_mesh_c[na] << '\n';
284 }
285 indn.n -= 2;
286 Iprintn(file, eesls->get_ees()->get_qe());
287 indn.n += 2;
288#ifdef USE_MEAN_COEF
289 Ifile << "ne, energy_mesh(ne), mean_coef_low_sigma:" << std::endl;
290 for (long ne = 0; ne < eesls->get_ees()->get_qe(); ne++) {
291 Ifile << std::setw(3) << ne << ' ' << std::setw(12)
292 << eesls->get_ees()->get_energy_mesh(ne) << " KeV " << std::setw(12)
293 << mean_coef_low_sigma[ne] << '\n';
294 }
295#else
296 Ifile << "ne, energy_mesh(ne), coef_low_sigma:" << std::endl;
297 for (long ne = 0; ne < eesls->get_ees()->get_qe(); ne++) {
298 Ifile << std::setw(3) << ne << ' ' << std::setw(12)
299 << eesls->get_ees()->get_energy_mesh(ne) << " KeV " << std::setw(12)
300 << coef_low_sigma[ne] << '\n';
301 }
302#endif
303 indn.n -= 2;
304}
305}
#define check_econd11(a, signb, stream)
Definition: FunNameStack.h:155
#define mfunname(string)
Definition: FunNameStack.h:45
PassivePtr< ElElasticScat > ees
virtual void print(std::ostream &file, int l) const
double get_sigma(double energy, double nscat) const
PassivePtr< HeedMatterDef > hmd
std::vector< PointsRan > low_angular_points_ran
HeedDeltaElectronCS()
Default constructor.
std::vector< double > beta
Table of velocities.
std::vector< double > momentum
Table of momenta [MeV/c].
std::vector< double > eLoss
std::vector< double > angular_mesh_c
Angular mesh, centers, angles in degrees.
std::vector< double > lambda
std::vector< PointsRan > angular_points_ran
PassivePtr< ElElasticScatLowSigma > eesls
static const double low_cut_angle_deg
std::vector< double > mean_coef_low_sigma
std::vector< double > low_lambda
static const long q_angular_mesh
Definition: BGMesh.cpp:5
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
double lorbeta(const double gamma_1)
as function of .
Definition: lorgamma.cpp:23
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:337
double e_cont_enloss(double ratio_Z_to_A, double I_eff, double density, double Ekin, double Ecut, double z)
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384
DoubleAc asin(const DoubleAc &f)
Definition: DoubleAc.cpp:470
indentation indn
Definition: prstream.cpp:15
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314
#define Ifile
Definition: prstream.h:196
#define mcerr
Definition: prstream.h:128
#define Iprintn(file, name)
Definition: prstream.h:205