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