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