Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
HeedMatterDef.cpp
Go to the documentation of this file.
1#include <fstream>
2#include <iomanip>
6
7// 2003, I. Smirnov
8
9namespace Heed {
10
11using CLHEP::pi;
12using CLHEP::mole;
13using CLHEP::gram;
14using CLHEP::cm3;
15using CLHEP::electron_mass_c2;
16using CLHEP::fine_structure_const;
17using CLHEP::Avogadro;
18
20 const std::vector<AtomPhotoAbsCS*>& faapacs,
21 double fW, double fF)
22 : matter(amatter), W(fW), F(fF), energy_mesh(fenergy_mesh) {
23 mfunname("HeedMatterDef::HeedMatterDef(...)");
24 check_econd11(matter, == nullptr, mcerr);
26 const long q = matter->qatom();
27 apacs.resize(q, nullptr);
28 for (long n = 0; n < q; ++n) {
29 apacs[n] = faapacs[n];
30 check_econd12(matter->atom(n)->Z(), !=, apacs[n]->get_Z(), mcerr);
31 }
32 check_econd11(F, == 0.0, mcerr);
33 if (W == 0.0) {
34#ifdef CALC_W_USING_CHARGES
35 double mean_I = 0.0;
36 double d = 0.0;
37 for (long n = 0; n < q; ++n) {
38 const double w = matter->weight_quan(n) * apacs[n]->get_Z();
39 mean_I += w * apacs[n]->get_I_min();
40 d += w;
41 }
42 W = coef_I_to_W * mean_I / d;
43#else
44 double mean_I = 0.0;
45 for (long n = 0; n < q; ++n) {
46 mean_I += matter->weight_quan(n) * apacs[n]->get_I_min();
47 }
48 W = coef_I_to_W * mean_I;
49#endif
50 }
51 inite_HeedMatterDef();
52}
53
55 const std::vector<MolecPhotoAbsCS*>& fampacs,
56 double fW, double fF)
57 : matter(agas), W(fW), F(fF), energy_mesh(fenergy_mesh) {
58 mfunname("HeedMatterDef::HeedMatterDef(...)");
59 check_econd11(agas, == nullptr, mcerr);
60 check_econd11(agas->qmolec(), <= 0, mcerr);
61 const long qat = agas->qatom();
62 apacs.resize(qat, nullptr);
63 const long qmol = agas->qmolec();
64 long nat = 0;
65 for (long nmol = 0; nmol < qmol; ++nmol) {
66 check_econd12(agas->molec(nmol)->tqatom(), !=, fampacs[nmol]->get_qatom(),
67 mcerr);
68 // number of different atoms in mol
69 const long qa = agas->molec(nmol)->qatom();
70 for (long na = 0; na < qa; ++na) {
71 apacs[nat] = fampacs[nmol]->get_atom(na);
72 check_econd12(apacs[nat]->get_Z(), !=, agas->molec(nmol)->atom(na)->Z(),
73 mcerr);
74 nat++;
75 }
76 }
77 if (F == 0.0) {
78#ifdef CALC_W_USING_CHARGES
79 double u = 0.0;
80 double d = 0.0;
81 for (long n = 0; n < qmol; ++n) {
82 const double w = agas->weight_quan_molec(n) * fampacs[n]->get_total_Z();
83 u += w * fampacs[n]->get_F();
84 d += w;
85 }
86 F = u / d;
87#else
88 for (long n = 0; n < qmol; ++n) {
89 F += agas->weight_quan_molec(n) * fampacs[n]->get_F();
90 }
91#endif
92 }
93
94 if (W == 0.0) {
95#ifdef CALC_W_USING_CHARGES
96 double u = 0.0;
97 double d = 0.0;
98 for (long n = 0; n < qmol; ++n) {
99 const double w = agas->weight_quan_molec(n) * fampacs[n]->get_total_Z();
100 u += w * fampacs[n]->get_W();
101 d += w;
102 }
103 W = u / d;
104#else
105 for (long n = 0; n < qmol; ++n) {
106 W += agas->weight_quan_molec(n) * fampacs[n]->get_W();
107 }
108#endif
109 }
110 inite_HeedMatterDef();
111}
112
114 const std::string& gas_notation,
115 const std::vector<MolecPhotoAbsCS*>& fampacs,
116 double fW, double fF)
117 : W(fW), F(fF), energy_mesh(fenergy_mesh) {
118 mfunnamep("HeedMatterDef::HeedMatterDef(...)");
119 MatterDef* amat = MatterDef::get_MatterDef(gas_notation);
120 GasDef* agas = dynamic_cast<GasDef*>(amat);
121 if (!agas) {
122 funnw.ehdr(mcerr);
123 mcerr << "notation supplied as the gas notation is not appear "
124 << "to be related to gas \n";
125 mcerr << "gas_notation=" << gas_notation << '\n';
126 spexit(mcerr);
127 }
128
129 matter = agas;
130 check_econd11(agas->qmolec(), <= 0, mcerr);
131 const long qat = agas->qatom();
132 apacs.resize(qat);
133 const long qmol = agas->qmolec();
134 long nat = 0;
135 for (long nmol = 0; nmol < qmol; ++nmol) {
136 check_econd12(agas->molec(nmol)->tqatom(), !=, fampacs[nmol]->get_qatom(),
137 mcerr);
138 // quantity of different atoms in molecule.
139 const long qa = agas->molec(nmol)->qatom();
140 for (long na = 0; na < qa; ++na) {
141 apacs[nat] = fampacs[nmol]->get_atom(na);
142 check_econd12(apacs[nat]->get_Z(), !=, agas->molec(nmol)->atom(na)->Z(),
143 mcerr);
144 nat++;
145 }
146 }
147 if (F == 0.0) {
148#ifdef CALC_W_USING_CHARGES
149 double u = 0.0;
150 double d = 0.0;
151 for (long n = 0; n < qmol; ++n) {
152 const double w = agas->weight_quan_molec(n) * fampacs[n]->get_total_Z();
153 u += w * fampacs[n]->get_F();
154 d += w;
155 }
156 F = u / d;
157#else
158 for (long n = 0; n < qmol; ++n) {
159 F += agas->weight_quan_molec(n) * fampacs[n]->get_F();
160 }
161#endif
162 }
163
164 if (W == 0.0) {
165#ifdef CALC_W_USING_CHARGES
166 double u = 0.0;
167 double d = 0.0;
168 for (long n = 0; n < qmol; ++n) {
169 const double w = agas->weight_quan_molec(n) * fampacs[n]->get_total_Z();
170 u += w * fampacs[n]->get_W();
171 d += w;
172 }
173 W = u / d;
174#else
175 for (long n = 0; n < qmol; ++n) {
176 W += agas->weight_quan_molec(n) * fampacs[n]->get_W();
177 }
178#endif
179 }
180 inite_HeedMatterDef();
181}
182
183void HeedMatterDef::inite_HeedMatterDef() {
184 mfunname("void HeedMatterDef::inite_HeedMatterDef()");
185 const double amean = matter->A_mean() / (gram / mole);
186 const double rho = matter->density() / (gram / cm3);
187 eldens_cm_3 = matter->Z_mean() / amean * Avogadro * rho;
190 wpla = eldens * 4. * pi * fine_structure_const / electron_mass_c2;
191 radiation_length = 0.0;
192 double rms = 0.0;
193 long qat = matter->qatom();
194 for (long n = 0; n < qat; ++n) {
195 rms += matter->atom(n)->A() * matter->weight_quan(n);
196 }
197 rms = rms / (gram / mole);
198
199 std::vector<double> RLenAt(qat);
200 std::vector<double> RuthAt(qat);
201 for (long n = 0; n < qat; ++n) {
202 const int z = matter->atom(n)->Z();
203 RLenAt[n] = 716.4 * matter->atom(n)->A() / (gram / mole) /
204 (z * (z + 1) * log(287. / sqrt(double(z))));
205 RuthAt[n] = 4. * pi * z * z * fine_structure_const * fine_structure_const;
206 }
207 std::vector<double> rm(qat);
208 for (long n = 0; n < qat; ++n) {
209 rm[n] = matter->atom(n)->A() / (gram / mole) * matter->weight_quan(n) / rms;
210 }
211 for (long n = 0; n < qat; ++n) {
212 radiation_length += rm[n] / RLenAt[n];
213 }
214 radiation_length = 1. / (rho * radiation_length);
215
216 Rutherford_const = 0.0;
217 for (long n = 0; n < qat; ++n) {
218 Rutherford_const += matter->weight_quan(n) * RuthAt[n];
219 }
220 Rutherford_const *= rho * Avogadro / amean;
221
222 min_ioniz_pot = DBL_MAX;
223 for (long n = 0; n < qat; ++n) {
224 if (min_ioniz_pot > apacs[n]->get_I_min()) {
225 min_ioniz_pot = apacs[n]->get_I_min();
226 }
227 }
228 long qe = energy_mesh->get_q();
229 ACS.resize(qe);
230 ICS.resize(qe);
231 epsip.resize(qe);
232 epsi1.resize(qe);
233 epsi2.resize(qe);
234 for (long ne = 0; ne < qe; ++ne) {
235 double e1 = energy_mesh->get_e(ne);
236 double e2 = energy_mesh->get_e(ne + 1);
237 double sa = 0.0;
238 double si = 0.0;
239 for (int na = 0; na < qat; ++na) {
240 const double ta = apacs[na]->get_integral_ACS(e1, e2);
241 sa += matter->weight_quan(na) * ta / (e2 - e1);
242 check_econd11a(ta, < 0, "ACS: ne=" << ne << " e1=" << e1 << " e2=" << e2
243 << " na=" << na << '\n',
244 mcerr);
245 const double ti =
247 ? apacs[na]->get_integral_TICS(e1, e2, min_ioniz_pot)
248 : apacs[na]->get_integral_ICS(e1, e2);
249 si += matter->weight_quan(na) * ti / (e2 - e1);
250 check_econd11a(ti, < 0, "ICS: ne=" << ne << " e1=" << e1 << " e2=" << e2
251 << " na=" << na << '\n',
252 mcerr);
253 }
254 ACS[ne] = sa;
255 check_econd11a(ACS[ne], < 0, "ne=" << ne << '\n', mcerr);
256 ICS[ne] = si;
257 check_econd11a(ICS[ne], < 0, "ne=" << ne << '\n', mcerr);
258
259 double ec = energy_mesh->get_ec(ne);
260 double ec2 = ec * ec;
261 epsip[ne] = -wpla / ec2;
262 epsi2[ne] = sa * C1_MEV2_MBN / ec * eldens / matter->Z_mean();
263 }
264
265 // To do next loop we need all epsi2
266 for (long ne = 0; ne < qe; ++ne) {
267 double ec = energy_mesh->get_ec(ne);
268 double ec2 = ec * ec;
269 double s = 0;
270 // integration of energy
271 for (long m = 0; m < qe; ++m) {
272 double em1 = energy_mesh->get_e(m);
273 double em2 = energy_mesh->get_e(m + 1);
274 double ecm = energy_mesh->get_ec(m);
275 if (m != ne) {
276 s += epsi2[m] * ecm * (em2 - em1) / (ecm * ecm - ec2);
277 } else {
278 double ee1 = 0.5 * (em1 + ecm);
279 double ee2 = 0.5 * (em2 + ecm);
280 double ep1, ep2; // extrapolated values to points ee1 and ee2
281 if (m == 0) {
282 ep1 = epsi2[m] + (ee1 - ecm) * (epsi2[m + 1] - epsi2[m]) /
283 (energy_mesh->get_ec(m + 1) - ecm);
284 } else {
285 ep1 = epsi2[m - 1] + (ee1 - energy_mesh->get_ec(m - 1)) *
286 (epsi2[m] - epsi2[m - 1]) /
287 (ecm - energy_mesh->get_ec(m - 1));
288 }
289 if (m < qe - 1) {
290 ep2 = epsi2[m] + (ee2 - ecm) * (epsi2[m + 1] - epsi2[m]) /
291 (energy_mesh->get_ec(m + 1) - ecm);
292 } else {
293 ep2 = epsi2[m] + (ee2 - ecm) * (epsi2[m] - epsi2[m - 1]) /
294 (ecm - energy_mesh->get_ec(m - 1));
295 }
296 s = s + ep1 * ee1 * (ecm - em1) / (ee1 * ee1 - ec2);
297 s = s + ep2 * ee2 * (em2 - ecm) / (ee2 * ee2 - ec2);
298 }
299 }
300 epsi1[ne] = (2. / pi) * s;
301 }
302}
303
304void HeedMatterDef::replace_epsi12(const std::string& file_name) {
305 mfunnamep("void HeedMatterDef::replace_epsi12(const std::string& file_name)");
306
307 std::ifstream file(file_name.c_str());
308 if (!file) {
309 funnw.ehdr(mcerr);
310 mcerr << "cannot open file " << file_name << std::endl;
311 spexit(mcerr);
312 } else {
313 mcout << "file " << file_name << " is opened" << std::endl;
314 }
315 long qe = 0; // number of points in input mesh
316 file >> qe;
317 check_econd11(qe, <= 2, mcerr);
318
319 std::vector<double> ener(qe);
320 std::vector<double> eps1(qe);
321 std::vector<double> eps2(qe);
322
323 for (long ne = 0; ne < qe; ++ne) {
324 file >> ener[ne] >> eps1[ne] >> eps2[ne];
325 check_econd11(eps2[ne], < 0.0, mcerr);
326 if (ne > 0) {
327 check_econd12(ener[ne], <, ener[ne - 1], mcerr);
328 }
329 }
330
332 double emin = ener[0] - 0.5 * (ener[1] - ener[0]);
333 double emax = ener[qe - 1] + 0.5 * (ener[qe - 1] - ener[qe - 2]);
334
335 qe = energy_mesh->get_q();
336 for (long ne = 0; ne < qe; ++ne) {
337 double ec = energy_mesh->get_ec(ne);
338 epsi1[ne] =
339 t_value_straight_point_ar<double, std::vector<double>,
341 pcmd, // dimension q
342 eps1, // dimension q-1
343 ec, 0, 1, emin, 1, emax);
344 epsi2[ne] =
345 t_value_straight_point_ar<double, std::vector<double>,
347 pcmd, // dimension q
348 eps2, // dimension q-1
349 ec, 1, 1, emin, 1, emax);
350 // Iprint3n(mcout, ec, epsi1[ne], epsi2[ne]);
351 }
352}
353
354void HeedMatterDef::print(std::ostream& file, int l) const {
355 if (l <= 0) return;
356 Ifile << "HeedMatterDef:\n";
357 indn.n += 2;
358 matter->print(file, 1);
359 if (l >= 2) {
360 long q = matter->qatom();
361 Ifile << "Printing " << q << " photoabsorption cross sections:\n";
362 indn.n += 2;
363 for (long n = 0; n < q; ++n) {
364 apacs[n]->print(file, l - 1);
365 }
366 indn.n -= 2;
367 }
368 Iprintan(file, eldens_cm_3, "1/cm^3");
369 Iprintan(file, eldens, "MeV^3");
370 Iprintan(file, xeldens, "MeV^2/cm");
371 Iprintn(file, wpla);
373 Iprintan(file, Rutherford_const, "1/cm^3");
374 Iprintn(file, W);
375 Iprintn(file, F);
376 Iprintn(file, min_ioniz_pot);
377 Iprintn(file, energy_mesh->get_q());
378 if (l >= 2) {
379 long qe = energy_mesh->get_q();
380 long ne;
381 indn.n += 2;
382 Ifile << " ne energy ACS(Mb) ICS(Mb) ACS(1/MeV^2) "
383 "ICS(1/MeV^2) epsip epsi1 epsi2 "
384 "(1+epsi1)^2+epsi2^2\n";
385 for (ne = 0; ne < qe; ne++) {
386 Ifile << std::setw(3) << ne << ' ' << std::setw(12)
387 << energy_mesh->get_e(ne) << ' ' << std::setw(12) << ACS[ne] << ' '
388 << std::setw(12) << ICS[ne] << ' ' << std::setw(12)
389 << ACS[ne] * C1_MEV2_MBN << ' ' << std::setw(12)
390 << ICS[ne] * C1_MEV2_MBN << ' ' << std::setw(12) << epsip[ne] << ' '
391 << std::setw(12) << epsi1[ne] << ' ' << std::setw(12) << epsi2[ne]
392 << ' ' << std::setw(12)
393 << pow((1 + epsi1[ne]), 2) + pow(epsi2[ne], 2) << " \n";
394 }
395 indn.n -= 2;
396 }
397 indn.n -= 2;
398}
399}
#define check_econd11(a, signb, stream)
Definition: FunNameStack.h:155
#define check_econd11a(a, signb, add, stream)
Definition: FunNameStack.h:172
#define mfunnamep(string)
Definition: FunNameStack.h:49
#define spexit(stream)
Definition: FunNameStack.h:256
#define check_econd12(a, sign, b, stream)
Definition: FunNameStack.h:163
#define mfunname(string)
Definition: FunNameStack.h:45
const std::vector< double > & weight_quan() const
Definition: AtomDef.h:136
long qatom() const
Definition: AtomDef.h:133
double Z_mean() const
Definition: AtomDef.h:140
const std::vector< AtomDef * > & atom() const
Definition: AtomDef.h:134
double A_mean() const
Definition: AtomDef.h:141
long get_q() const
Return number of bins.
Definition: EnergyMesh.h:42
double get_ec(long n) const
Return center of a given bin.
Definition: EnergyMesh.h:50
double get_e(long n) const
Return left side of a given bin.
Definition: EnergyMesh.h:48
long qmolec() const
Definition: GasDef.h:47
const std::vector< double > & weight_quan_molec() const
Definition: GasDef.h:52
const std::vector< MoleculeDef * > & molec() const
Definition: GasDef.h:48
double eldens
Electron density MeV**3.
Definition: HeedMatterDef.h:32
double W
Mean work per pair production, MeV.
Definition: HeedMatterDef.h:37
double wpla
Squared plasma energy;.
Definition: HeedMatterDef.h:34
double xeldens
Long. electron density MeV**2/cm (for x=1 cm).
Definition: HeedMatterDef.h:33
std::vector< const AtomPhotoAbsCS * > apacs
Definition: HeedMatterDef.h:29
HeedMatterDef()=default
Default constructor.
double eldens_cm_3
Electron density cm**-3.
Definition: HeedMatterDef.h:31
std::vector< double > ACS
Photoabsorption cross section per one atom(Mb).
Definition: HeedMatterDef.h:43
std::vector< double > epsip
Some plasma dielectric constant (not used, but just initialized for print)
Definition: HeedMatterDef.h:48
MatterDef * matter
Definition: HeedMatterDef.h:28
std::vector< double > epsi1
Real part of dielectric constant (e_1 - 1).
Definition: HeedMatterDef.h:50
double Rutherford_const
Const for Rutherford cross section (1/cm3).
Definition: HeedMatterDef.h:36
std::vector< double > ICS
Definition: HeedMatterDef.h:44
std::vector< double > epsi2
Imaginary part of dielectric constant.
Definition: HeedMatterDef.h:52
double F
Fano factor.
Definition: HeedMatterDef.h:38
double radiation_length
Radiation Length.
Definition: HeedMatterDef.h:35
void print(std::ostream &file, int l) const
EnergyMesh * energy_mesh
Definition: HeedMatterDef.h:39
void replace_epsi12(const std::string &file_name)
static constexpr int s_use_mixture_thresholds
Definition: HeedMatterDef.h:90
static MatterDef * get_MatterDef(const std::string &fnotation)
Definition: MatterDef.cpp:120
void print(std::ostream &file, int l) const
Definition: MatterDef.cpp:100
double density() const
Definition: MatterDef.h:51
Definition: BGMesh.cpp:6
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:337
constexpr double C1_MEV2_MBN
constexpr double coef_I_to_W
Definition: PhotoAbsCS.h:585
constexpr double C1_MEV_CM
indentation indn
Definition: prstream.cpp:15
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314
#define Iprintan(file, name, addition)
Definition: prstream.h:211
#define mcout
Definition: prstream.h:126
#define Ifile
Definition: prstream.h:195
#define mcerr
Definition: prstream.h:128
#define Iprintn(file, name)
Definition: prstream.h:204