Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
MediumGas.cc
Go to the documentation of this file.
1#include <algorithm>
2#include <cstdio>
3#include <cstdlib>
4#include <cstring>
5#include <ctime>
6#include <fstream>
7#include <iomanip>
8#include <iostream>
9#include <sstream>
10#include <numeric>
11#include <utility>
12
15#include "Garfield/MediumGas.hh"
17#include "Garfield/Utilities.hh"
18
19namespace {
20
21std::string FmtFloat(const double x, const unsigned int width = 15,
22 const unsigned int precision = 8) {
23 char buffer[256];
24 std::snprintf(buffer, width + 1, "%*.*E", width, precision, x);
25 return std::string(buffer);
26}
27
28std::string FmtInt(const int n, const unsigned int width) {
29 char buffer[256];
30 std::snprintf(buffer, width + 1, "%*d", width, n);
31 return std::string(buffer);
32}
33
34void PrintArray(const std::vector<double>& values, std::ofstream& outfile,
35 int& col, const int ncols) {
36 for (const auto value : values) {
37 outfile << FmtFloat(value);
38 ++col;
39 if (col % ncols == 0) outfile << "\n";
40 }
41}
42
43void PrintExtrapolation(const std::pair<unsigned int, unsigned int>& extr) {
44 std::cout << " Low field extrapolation: ";
45 if (extr.first == 0)
46 std::cout << " constant\n";
47 else if (extr.first == 1)
48 std::cout << " linear\n";
49 else if (extr.first == 2)
50 std::cout << " exponential\n";
51 else
52 std::cout << " unknown\n";
53 std::cout << " High field extrapolation: ";
54 if (extr.second == 0)
55 std::cout << " constant\n";
56 else if (extr.second == 1)
57 std::cout << " linear\n";
58 else if (extr.second == 2)
59 std::cout << " exponential\n";
60 else
61 std::cout << " unknown\n";
62}
63
64void PrintAbsentInNew(const std::string& par) {
65 std::cout << " Warning: The " << par << " is absent in the dataset "
66 << "to be added; data reset.\n";
67}
68
69void PrintAbsentInExisting(const std::string& par) {
70 std::cout << " Warning: The " << par << " is absent in the existing data; "
71 << "new data not used.\n";
72}
73
74bool Similar(const double v1, const double v2, const double eps) {
75 const double dif = v1 - v2;
76 const double sum = fabs(v1) + fabs(v2);
77 return fabs(dif) < std::max(eps * sum, Garfield::Small);
78}
79
80int Equal(const std::vector<double>& fields1,
81 const std::vector<double>& fields2, const double eps) {
82 if (fields1.size() != fields2.size()) return 0;
83 const size_t n = fields1.size();
84 for (size_t i = 0; i < n; ++i) {
85 if (!Similar(fields1[i], fields2[i], eps)) return 0;
86 }
87 return 1;
88}
89
90int FindIndex(const double field, const std::vector<double>& fields,
91 const double eps) {
92
93 const int n = fields.size();
94 for (int i = 0; i < n; ++i) {
95 if (Similar(field, fields[i], eps)) return i;
96 }
97 return -1;
98}
99
100void ResizeA(std::vector<std::vector<std::vector<double> > >& tab,
101 const int ne, const int nb, const int na) {
102 if (tab.empty()) return;
103 tab.resize(na,
104 std::vector<std::vector<double> >(nb,
105 std::vector<double>(ne, 0.)));
106}
107
108}
109
110namespace Garfield {
111
113 : Medium(), m_pressureTable(m_pressure), m_temperatureTable(m_temperature) {
114 m_className = "MediumGas";
115
116 m_gas.fill("");
117 m_fraction.fill(0.);
118 m_atWeight.fill(0.);
119 m_atNum.fill(0.);
120 // Default gas mixture: pure argon
121 m_gas[0] = "Ar";
122 m_fraction[0] = 1.;
123 m_name = m_gas[0];
125
126 m_rPenningGas.fill(0.);
127 m_lambdaPenningGas.fill(0.);
128
129 m_isChanged = true;
130
131 EnableDrift();
133}
134
135bool MediumGas::SetComposition(const std::string& gas1, const double f1,
136 const std::string& gas2, const double f2,
137 const std::string& gas3, const double f3,
138 const std::string& gas4, const double f4,
139 const std::string& gas5, const double f5,
140 const std::string& gas6, const double f6) {
141 std::array<std::string, 6> gases = {gas1, gas2, gas3, gas4, gas5, gas6};
142 std::array<double, 6> fractions = {f1, f2, f3, f4, f5, f6};
143
144 // Make a backup copy of the gas composition.
145 const std::array<std::string, m_nMaxGases> gasOld = m_gas;
146 const unsigned int nComponentsOld = m_nComponents;
147
148 // Reset all transport parameters.
149 ResetTables();
150 // Force a recalculation of the collision rates.
151 m_isChanged = true;
152
153 m_nComponents = 0;
154 m_gas.fill("");
155 m_fraction.fill(0.);
156 m_atWeight.fill(0.);
157 m_atNum.fill(0.);
158 for (unsigned int i = 0; i < 6; ++i) {
159 if (fractions[i] < Small) continue;
160 // Find the gas name corresponding to the input string.
161 const std::string gasname = GetGasName(gases[i]);
162 if (!gasname.empty()) {
163 m_gas[m_nComponents] = gasname;
164 m_fraction[m_nComponents] = fractions[i];
166 }
167 }
168
169 // Check if at least one valid ingredient was specified.
170 if (m_nComponents == 0) {
171 std::cerr << m_className << "::SetComposition:\n"
172 << " Error setting the composition. No valid components.\n";
173 return false;
174 }
175
176 // Establish the name.
177 m_name = "";
178 double sum = 0.;
179 for (unsigned int i = 0; i < m_nComponents; ++i) {
180 if (i > 0) m_name += "/";
181 m_name += m_gas[i];
182 sum += m_fraction[i];
183 }
184 // Normalise the fractions to one.
185 for (unsigned int i = 0; i < m_nComponents; ++i) {
186 m_fraction[i] /= sum;
187 }
188
189 // Set the atomic weight and number.
190 for (unsigned int i = 0; i < m_nComponents; ++i) {
192 }
193
194 // Print the composition.
195 std::cout << m_className << "::SetComposition:\n " << m_name;
196 if (m_nComponents > 1) {
197 std::cout << " (" << m_fraction[0] * 100;
198 for (unsigned int i = 1; i < m_nComponents; ++i) {
199 std::cout << "/" << m_fraction[i] * 100;
200 }
201 std::cout << ")";
202 }
203 std::cout << "\n";
204
205
206 // Copy the previous Penning transfer parameters.
207 std::array<double, m_nMaxGases> rPenningGasOld;
208 std::array<double, m_nMaxGases> lambdaPenningGasOld;
209 rPenningGasOld.fill(0.);
210 lambdaPenningGasOld.fill(0.);
211 rPenningGasOld.swap(m_rPenningGas);
212 lambdaPenningGasOld.swap(m_lambdaPenningGas);
213 for (unsigned int i = 0; i < m_nComponents; ++i) {
214 for (unsigned int j = 0; j < nComponentsOld; ++j) {
215 if (m_gas[i] != gasOld[j]) continue;
216 if (rPenningGasOld[j] < Small) continue;
217 m_rPenningGas[i] = rPenningGasOld[j];
218 m_lambdaPenningGas[i] = lambdaPenningGasOld[i];
219 std::cout << m_className << "::SetComposition:\n"
220 << " Using Penning transfer parameters for " << m_gas[i]
221 << " from previous mixture.\n"
222 << " r = " << m_rPenningGas[i] << "\n"
223 << " lambda = " << m_lambdaPenningGas[i] << " cm\n";
224 }
225 }
226 return true;
227}
228
229void MediumGas::GetComposition(std::string& gas1, double& f1, std::string& gas2,
230 double& f2, std::string& gas3, double& f3,
231 std::string& gas4, double& f4, std::string& gas5,
232 double& f5, std::string& gas6, double& f6) {
233 gas1 = m_gas[0];
234 gas2 = m_gas[1];
235 gas3 = m_gas[2];
236 gas4 = m_gas[3];
237 gas5 = m_gas[4];
238 gas6 = m_gas[5];
239 f1 = m_fraction[0];
240 f2 = m_fraction[1];
241 f3 = m_fraction[2];
242 f4 = m_fraction[3];
243 f5 = m_fraction[4];
244 f6 = m_fraction[5];
245}
246
247void MediumGas::GetComponent(const unsigned int i, std::string& label,
248 double& f) {
249 if (i >= m_nComponents) {
250 std::cerr << m_className << "::GetComponent: Index out of range.\n";
251 label = "";
252 f = 0.;
253 return;
254 }
255
256 label = m_gas[i];
257 f = m_fraction[i];
258}
259
260void MediumGas::SetAtomicNumber(const double /*z*/) {
261 std::cerr << m_className << "::SetAtomicNumber:\n"
262 << " Effective Z cannot be changed directly.\n"
263 << " Use SetComposition to define the gas mixture.\n";
264}
265
266void MediumGas::SetAtomicWeight(const double /*a*/) {
267 std::cerr << m_className << "::SetAtomicWeight:\n"
268 << " Effective A cannot be changed directly.\n"
269 << " Use SetComposition to define the gas mixture.\n";
270}
271
272void MediumGas::SetNumberDensity(const double /*n*/) {
273 std::cerr << m_className << "::SetNumberDensity:\n"
274 << " Density cannot be changed directly.\n"
275 << " Use SetTemperature and SetPressure.\n";
276}
277
278void MediumGas::SetMassDensity(const double /*rho*/) {
279 std::cerr << m_className << "::SetMassDensity:\n"
280 << " Density cannot be changed directly.\n"
281 << " Use SetTemperature, SetPressure and SetComposition.\n";
282}
283
285 // Effective A, weighted by the fractions of the components.
286 double a = 0.;
287 for (unsigned int i = 0; i < m_nComponents; ++i) {
288 a += m_atWeight[i] * m_fraction[i];
289 }
290 return a;
291}
292
294 // Ideal gas law.
295 return LoschmidtNumber * (m_pressure / AtmosphericPressure) *
296 (ZeroCelsius / m_temperature);
297}
298
300 return GetNumberDensity() * GetAtomicWeight() * AtomicMassUnit;
301}
302
304 // Effective Z, weighted by the fractions of the components.
305 double z = 0.;
306 for (unsigned int i = 0; i < m_nComponents; ++i) {
307 z += m_atNum[i] * m_fraction[i];
308 }
309 return z;
310}
311
312bool MediumGas::LoadGasFile(const std::string& filename) {
313
314 // -----------------------------------------------------------------------
315 // GASGET
316 // -----------------------------------------------------------------------
317
318 std::ifstream gasfile;
319 // Open the file.
320 gasfile.open(filename.c_str());
321 // Make sure the file could be opened.
322 if (!gasfile.is_open()) {
323 std::cerr << m_className << "::LoadGasFile:\n"
324 << " Cannot open file " << filename << ".\n";
325 return false;
326 }
327 std::cout << m_className << "::LoadGasFile: Reading " << filename << ".\n";
328
329 ResetTables();
330
331 // Start reading the data.
332 if (m_debug) std::cout << m_className << "::LoadGasFile: Reading header.\n";
333 int version = 12;
334 // GASOK bits
335 std::bitset<20> gasok;
336 // Gas composition
337 constexpr int nMagboltzGases = 60;
338 std::vector<double> mixture(nMagboltzGases, 0.);
339 if (!ReadHeader(gasfile, version, gasok, m_tab2d, mixture,
341 gasfile.close();
342 return false;
343 }
344 std::cout << m_className << "::LoadGasFile: Version " << version << "\n";
345
346 // Check the gas mixture.
347 std::vector<std::string> gasnames;
348 std::vector<double> percentages;
349 if (!GetMixture(mixture, version, gasnames, percentages)) {
350 std::cerr << m_className << "::LoadGasFile:\n "
351 << "Cannot determine the gas composition.\n";
352 gasfile.close();
353 return false;
354 }
355
356 m_name = "";
357 m_nComponents = gasnames.size();
358 for (unsigned int i = 0; i < m_nComponents; ++i) {
359 if (i > 0) m_name += "/";
360 m_name += gasnames[i];
361 m_gas[i] = gasnames[i];
362 m_fraction[i] = percentages[i] / 100.;
364 }
365 std::cout << m_className << "::LoadGasFile:\n"
366 << " Gas composition set to " << m_name;
367 if (m_nComponents > 1) {
368 std::cout << " (" << m_fraction[0] * 100;
369 for (unsigned int i = 1; i < m_nComponents; ++i) {
370 std::cout << "/" << m_fraction[i] * 100;
371 }
372 std::cout << ")";
373 }
374 std::cout << "\n";
375
376 const int nE = m_eFields.size();
377 const int nB = m_bFields.size();
378 const int nA = m_bAngles.size();
379 if (m_debug) {
380 std::cout << m_className << "::LoadGasFile:\n " << nE
381 << " electric field(s), " << nB
382 << " magnetic field(s), " << nA << " angle(s).\n";
383 }
384
385 // Decode the GASOK bits.
386 // GASOK(I) : .TRUE. if present
387 // (1) electron drift velocity || E
388 // (2) ion mobility,
389 // (3) longitudinal diffusion || E
390 // (4) Townsend coefficient,
391 // (5) cluster size distribution.
392 // (6) attachment coefficient,
393 // (7) Lorentz angle,
394 // (8) transverse diffusion || ExB and Bt
395 // (9) electron drift velocity || Bt
396 // (10) electron drift velocity || ExB
397 // (11) diffusion tensor
398 // (12) ion dissociation
399 // (13) allocated for SRIM data (not used)
400 // (14) allocated for HEED data (not used)
401 // (15) excitation rates
402 // (16) ionisation rates
403
404 if (gasok[0]) Init(nE, nB, nA, m_eVelE, 0.);
405 if (gasok[1]) Init(nE, nB, nA, m_iMob, 0.);
406 if (gasok[2]) Init(nE, nB, nA, m_eDifL, 0.);
407 if (gasok[3]) {
408 Init(nE, nB, nA, m_eAlp, -30.);
409 Init(nE, nB, nA, m_eAlp0, -30.);
410 }
411 if (gasok[5]) Init(nE, nB, nA, m_eAtt, -30.);
412 if (gasok[6]) Init(nE, nB, nA, m_eLor, -30.);
413 if (gasok[7]) Init(nE, nB, nA, m_eDifT, 0.);
414 if (gasok[8]) Init(nE, nB, nA, m_eVelB, 0.);
415 if (gasok[9]) Init(nE, nB, nA, m_eVelX, 0.);
416 if (gasok[10]) Init(nE, nB, nA, 6, m_eDifM, 0.);
417 if (gasok[11]) Init(nE, nB, nA, m_iDis, -30.);
418 if (gasok[14]) Init(nE, nB, nA, m_excLevels.size(), m_excRates, 0.);
419 if (gasok[15]) Init(nE, nB, nA, m_ionLevels.size(), m_ionRates, 0.);
420
421 // Force re-initialisation of collision rates etc.
422 m_isChanged = true;
423
424 if (m_debug) {
425 const std::string fmt = m_tab2d ? "3D" : "1D";
426 std::cout << m_className << "::LoadGasFile: Reading " << fmt << " table.\n";
427 }
428
429 // Drift velocity along E, Bt and ExB
430 double ve = 0., vb = 0., vx = 0.;
431 // Lorentz angle
432 double lor = 0.;
433 // Longitudinal and transverse diffusion coefficients
434 double dl = 0., dt = 0.;
435 // Townsend and attachment coefficients
436 double alpha = 0., alpha0 = 0., eta = 0.;
437 // Ion mobility and dissociation coefficient
438 double mu = 0., dis = 0.;
439 // Diffusion tensor.
440 std::array<double, 6> diff;
441 // Excitation and ionization rates.
442 const unsigned int nexc = m_excLevels.size();
443 std::vector<double> rexc(nexc, 0.);
444 const unsigned int nion = m_ionLevels.size();
445 std::vector<double> rion(nion, 0.);
446 for (int i = 0; i < nE; i++) {
447 for (int j = 0; j < nA; j++) {
448 for (int k = 0; k < nB; k++) {
449 if (m_tab2d) {
450 ReadRecord3D(gasfile, ve, vb, vx, dl, dt, alpha, alpha0, eta, mu,
451 lor, dis, diff, rexc, rion);
452 } else {
453 ReadRecord1D(gasfile, ve, vb, vx, dl, dt, alpha, alpha0, eta, mu,
454 lor, dis, diff, rexc, rion);
455 }
456 if (!m_eVelE.empty()) m_eVelE[j][k][i] = ve;
457 if (!m_eVelB.empty()) m_eVelB[j][k][i] = vb;
458 if (!m_eVelX.empty()) m_eVelX[j][k][i] = vx;
459 if (!m_eDifL.empty()) m_eDifL[j][k][i] = dl;
460 if (!m_eDifT.empty()) m_eDifT[j][k][i] = dt;
461 if (!m_eAlp.empty()) {
462 m_eAlp[j][k][i] = alpha;
463 m_eAlp0[j][k][i] = alpha0;
464 }
465 if (!m_eAtt.empty()) m_eAtt[j][k][i] = eta;
466 if (!m_iMob.empty()) m_iMob[j][k][i] = mu;
467 if (!m_eLor.empty()) m_eLor[j][k][i] = lor;
468 if (!m_iDis.empty()) m_iDis[j][k][i] = dis;
469 if (!m_eDifM.empty()) {
470 for (int l = 0; l < 6; l++) {
471 m_eDifM[l][j][k][i] = diff[l];
472 }
473 }
474 if (!m_excRates.empty()) {
475 for (unsigned int l = 0; l < nexc; ++l) {
476 m_excRates[l][j][k][i] = rexc[l];
477 }
478 }
479 if (!m_ionRates.empty()) {
480 for (unsigned int l = 0; l < nion; ++l) {
481 m_ionRates[l][j][k][i] = rion[l];
482 }
483 }
484 }
485 }
486 }
487
488 // Extrapolation methods
489 std::array<unsigned int, 13> extrapH = {{0}};
490 std::array<unsigned int, 13> extrapL = {{1}};
491 // Interpolation methods
492 std::array<unsigned int, 13> interp = {{2}};
493 // Ion diffusion coefficients.
494 double ionDiffL = 0.;
495 double ionDiffT = 0.;
496 // Gas pressure [Torr] and temperature [K].
497 double pgas = 0.;
498 double tgas = 0.;
499 // Moving on to the file footer
500 gasfile.ignore(std::numeric_limits<std::streamsize>::max(), '\n');
501 if (m_debug) std::cout << m_className << "::LoadGasFile: Reading footer.\n";
502 ReadFooter(gasfile, extrapH, extrapL, interp,
503 m_eThrAlp, m_eThrAtt, m_iThrDis, ionDiffL, ionDiffT, pgas, tgas);
504 gasfile.close();
505
506 // Decrement the threshold indices (compatibility with Fortran).
507 if (m_eThrAlp > 0) --m_eThrAlp;
508 if (m_eThrAtt > 0) --m_eThrAtt;
509 if (m_iThrDis > 0) --m_iThrDis;
510
511 // Set the reference pressure and temperature.
512 if (pgas > 0.) m_pressure = m_pressureTable = pgas;
513 if (tgas > 0.) m_temperature = m_temperatureTable = tgas;
514
515 // Multiply the E/p values by the pressure.
516 for (auto& field : m_eFields) field *= m_pressureTable;
517
518 // Scale the parameters.
519 const double sqrp = sqrt(m_pressureTable);
520 const double logp = log(m_pressureTable);
521 for (int i = nE; i--;) {
522 for (int j = nA; j--;) {
523 for (int k = nB; k--;) {
524 if (!m_eDifL.empty()) m_eDifL[j][k][i] /= sqrp;
525 if (!m_eDifT.empty()) m_eDifT[j][k][i] /= sqrp;
526 if (!m_eDifM.empty()) {
527 for (int l = 6; l--;) m_eDifM[l][j][k][i] /= m_pressureTable;
528 }
529 if (!m_eAlp.empty()) {
530 m_eAlp[j][k][i] += logp;
531 m_eAlp0[j][k][i] += logp;
532 }
533 if (!m_eAtt.empty()) m_eAtt[j][k][i] += logp;
534 if (!m_iDis.empty()) m_iDis[j][k][i] += logp;
535 /*
536 for (auto& exc : m_excRates) {
537 exc[j][k][i] /= m_pressureTable;
538 }
539 for (auto& ion : m_ionRates) {
540 ion[j][k][i] /= m_pressureTable;
541 }
542 */
543 }
544 }
545 }
546
547 // Decode the extrapolation and interpolation tables.
548 m_extrVel = {extrapL[0], extrapH[0]};
549 // Indices 1 and 2 correspond to velocities along Bt and ExB.
550 m_extrDif = {extrapL[3], extrapH[3]};
551 m_extrAlp = {extrapL[4], extrapH[4]};
552 m_extrAtt = {extrapL[5], extrapH[5]};
553 m_extrMob = {extrapL[6], extrapH[6]};
554 m_extrLor = {extrapL[7], extrapH[7]};
555 // Index 8: transverse diffusion.
556 m_extrDis = {extrapL[9], extrapH[9]};
557 // Index 10: diff. tensor
558 m_extrExc = {extrapL[11], extrapH[11]};
559 m_extrIon = {extrapL[12], extrapH[12]};
560 m_intpVel = interp[0];
561 m_intpDif = interp[3];
562 m_intpAlp = interp[4];
563 m_intpAtt = interp[5];
564 m_intpMob = interp[6];
565 m_intpLor = interp[7];
566 m_intpDis = interp[9];
567 m_intpExc = interp[11];
568 m_intpIon = interp[12];
569
570 // Ion diffusion
571 if (ionDiffL > 0.) Init(nE, nB, nA, m_iDifL, ionDiffL);
572 if (ionDiffT > 0.) Init(nE, nB, nA, m_iDifT, ionDiffT);
573
574 if (m_debug) std::cout << m_className << "::LoadGasFile: Done.\n";
575
576 return true;
577}
578
579bool MediumGas::ReadHeader(std::ifstream& gasfile, int& version,
580 std::bitset<20>& gasok, bool& is3d, std::vector<double>& mixture,
581 std::vector<double>& efields, std::vector<double>& bfields,
582 std::vector<double>& angles, std::vector<ExcLevel>& excLevels,
583 std::vector<IonLevel>& ionLevels) {
584
585 gasok.reset();
586 bool done = false;
587 char line[256];
588 while (gasfile.getline(line, 256)) {
589 const bool quotes = (strstr(line, "\"") != NULL);
590 if (strncmp(line, " The gas tables follow:", 8) == 0 ||
591 strncmp(line, "The gas tables follow:", 7) == 0) {
592 done = true;
593 break;
594 }
595 char* token = strtok(line, " :,%");
596 while (token) {
597 if (strcmp(token, "Version") == 0) {
598 token = strtok(NULL, " :,%");
599 version = atoi(token);
600 // Check the version number.
601 if (version != 10 && version != 11 && version != 12) {
602 std::cerr << m_className << "::ReadHeader:\n"
603 << " The file has version number " << version << ".\n"
604 << " Files written in this format cannot be read.\n";
605 return false;
606 }
607 } else if (strcmp(token, "GASOK") == 0) {
608 // Get the GASOK bits indicating if a parameter
609 // is present in the table (T) or not (F).
610 token = strtok(NULL, " :,%\t");
611 token = strtok(NULL, " :,%\t");
612 std::string okstr(token);
613 if (m_debug) std::cout << " GASOK bits: " << okstr << "\n";
614 if (okstr.size() < 20) {
615 std::cerr << m_className << "::ReadHeader:\n"
616 << " Unexpected size of GASOK string ("
617 << okstr.size() << ").\n";
618 return false;
619 }
620 for (unsigned int i = 0; i < 20; ++i) {
621 if (okstr[i] == 'T') gasok.set(i);
622 }
623 } else if (strcmp(token, "Identifier") == 0) {
624 // Get the identification string.
625 std::string identifier = "";
626 token = strtok(NULL, "\n");
627 if (token != NULL) identifier += token;
628 if (m_debug) std::cout << " Identifier: " << identifier << "\n";
629 } else if (strcmp(token, "Dimension") == 0) {
630 token = strtok(NULL, " :,%\t");
631 if (strcmp(token, "F") == 0) {
632 is3d = false;
633 } else {
634 is3d = true;
635 }
636 token = strtok(NULL, " :,%\t");
637 const int nE = atoi(token);
638 // Check the number of E points.
639 if (nE <= 0) {
640 std::cerr << m_className << "::ReadHeader:\n"
641 << " Number of E fields out of range.\n";
642 return false;
643 }
644 token = strtok(NULL, " :,%\t");
645 const int nA = atoi(token);
646 // Check the number of angles.
647 if (is3d && nA <= 0) {
648 std::cerr << m_className << "::ReadHeader:\n"
649 << " Number of E-B angles out of range.\n";
650 return false;
651 }
652 token = strtok(NULL, " :,%\t");
653 const int nB = atoi(token);
654 // Check the number of B points.
655 if (is3d && nB <= 0) {
656 std::cerr << m_className << "::ReadHeader:\n"
657 << " Number of B fields out of range.\n";
658 return false;
659 }
660 efields.resize(nE);
661 angles.resize(nA);
662 bfields.resize(nB);
663 // Fill in the excitation/ionisation structs
664 // Excitation
665 token = strtok(NULL, " :,%\t");
666 const int nexc = atoi(token);
667 // Ionization
668 token = strtok(NULL, " :,%\t");
669 const int nion = atoi(token);
670 if (m_debug) {
671 std::cout << " " << nexc << " excitations, " << nion
672 << " ionisations.\n";
673 }
674 } else if (strcmp(token, "E") == 0) {
675 token = strtok(NULL, " :,%");
676 if (strncmp(token, "fields", 6) == 0) {
677 const int nE = efields.size();
678 for (int i = 0; i < nE; ++i) gasfile >> efields[i];
679 }
680 } else if (strcmp(token, "E-B") == 0) {
681 token = strtok(NULL, " :,%");
682 if (strncmp(token, "angles", 6) == 0) {
683 const int nA = angles.size();
684 for (int i = 0; i < nA; ++i) gasfile >> angles[i];
685 }
686 } else if (strcmp(token, "B") == 0) {
687 token = strtok(NULL, " :,%");
688 if (strncmp(token, "fields", 6) == 0) {
689 double bstore = 0.;
690 const int nB = bfields.size();
691 for (int i = 0; i < nB; i++) {
692 gasfile >> bstore;
693 bfields[i] = bstore / 100.;
694 }
695 }
696 } else if (strcmp(token, "Mixture") == 0) {
697 const unsigned int nMagboltzGases = mixture.size();
698 for (unsigned int i = 0; i < nMagboltzGases; ++i) {
699 gasfile >> mixture[i];
700 }
701 } else if (strcmp(token, "Excitation") == 0) {
702 // Skip number.
703 token = strtok(NULL, " :,%");
704 ExcLevel exc;
705 // Get label.
706 if (quotes) {
707 token = strtok(NULL, "\"");
708 token = strtok(NULL, "\"");
709 } else {
710 token = strtok(NULL, " :,%");
711 }
712 exc.label = token;
713 // Get energy.
714 token = strtok(NULL, " :,%");
715 exc.energy = atof(token);
716 // Get Penning probability.
717 token = strtok(NULL, " :,%");
718 exc.prob = atof(token);
719 exc.rms = 0.;
720 exc.dt = 0.;
721 if (version >= 11) {
722 // Get Penning rms distance.
723 token = strtok(NULL, " :,%");
724 if (token) {
725 exc.rms = atof(token);
726 // Get decay time.
727 token = strtok(NULL, " :,%");
728 if (token) exc.dt = atof(token);
729 }
730 }
731 excLevels.push_back(std::move(exc));
732 } else if (strcmp(token, "Ionisation") == 0) {
733 // Skip number.
734 token = strtok(NULL, " :,%");
735 // Get label.
736 IonLevel ion;
737 if (quotes) {
738 token = strtok(NULL, "\"");
739 token = strtok(NULL, "\"");
740 } else {
741 token = strtok(NULL, " :,%");
742 }
743 ion.label += token;
744 // Get energy.
745 token = strtok(NULL, " :,%");
746 ion.energy = atof(token);
747 ionLevels.push_back(std::move(ion));
748 }
749 token = strtok(NULL, " :,%");
750 }
751 }
752 return done;
753}
754
755void MediumGas::ReadRecord3D(std::ifstream& gasfile,
756 double& ve, double& vb, double& vx, double& dl, double& dt,
757 double& alpha, double& alpha0, double& eta, double& mu, double& lor,
758 double& dis, std::array<double, 6>& dif,
759 std::vector<double>& rexc, std::vector<double>& rion) {
760
761 // Drift velocity along E, Bt and ExB
762 gasfile >> ve >> vb >> vx;
763 // Convert from cm / us to cm / ns.
764 ve *= 1.e-3;
765 vb *= 1.e-3;
766 vx *= 1.e-3;
767 // Longitudinal and transverse diffusion coefficients
768 gasfile >> dl >> dt;
769 // Townsend and attachment coefficients
770 gasfile >> alpha >> alpha0 >> eta;
771 // Ion mobility
772 gasfile >> mu;
773 // Convert from cm2 / (V us) to cm2 / (V ns)
774 mu *= 1.e-3;
775 // Lorentz angle
776 gasfile >> lor;
777 // Ion dissociation
778 gasfile >> dis;
779 // Diffusion tensor
780 for (int l = 0; l < 6; l++) gasfile >> dif[l];
781 // Excitation rates
782 const unsigned int nexc = rexc.size();
783 for (unsigned int l = 0; l < nexc; ++l) gasfile >> rexc[l];
784 // Ionization rates
785 const unsigned int nion = rion.size();
786 for (unsigned int l = 0; l < nion; ++l) gasfile >> rion[l];
787}
788
789void MediumGas::ReadRecord1D(std::ifstream& gasfile,
790 double& ve, double& vb, double& vx, double& dl, double& dt,
791 double& alpha, double& alpha0, double& eta, double& mu, double& lor,
792 double& dis, std::array<double, 6>& dif,
793 std::vector<double>& rexc, std::vector<double>& rion) {
794
795 double waste = 0.;
796 gasfile >> ve >> waste >> vb >> waste >> vx >> waste;
797 ve *= 1.e-3;
798 vb *= 1.e-3;
799 vx *= 1.e-3;
800 gasfile >> dl >> waste >> dt >> waste;
801 gasfile >> alpha >> waste >> alpha0 >> eta >> waste;
802 gasfile >> mu >> waste;
803 mu *= 1.e-3;
804 gasfile >> lor >> waste;
805 gasfile >> dis >> waste;
806 for (int j = 0; j < 6; j++) gasfile >> dif[j] >> waste;
807 const unsigned int nexc = rexc.size();
808 for (unsigned int j = 0; j < nexc; ++j) gasfile >> rexc[j] >> waste;
809 const unsigned int nion = rion.size();
810 for (unsigned int j = 0; j < nion; ++j) gasfile >> rion[j] >> waste;
811}
812
813void MediumGas::ReadFooter(std::ifstream& gasfile,
814 std::array<unsigned int, 13>& extrapH,
815 std::array<unsigned int, 13>& extrapL,
816 std::array<unsigned int, 13>& interp,
817 unsigned int& thrAlp, unsigned int& thrAtt, unsigned int& thrDis,
818 double& ionDiffL, double& ionDiffT,
819 double& pgas, double& tgas) {
820
821 bool done = false;
822 while (!done) {
823 char line[256];
824 gasfile.getline(line, 256);
825 char* token = strtok(line, " :,%=\t");
826 while (token) {
827 if (strcmp(token, "H") == 0) {
828 token = strtok(NULL, " :,%=\t");
829 for (int i = 0; i < 13; i++) {
830 token = strtok(NULL, " :,%=\t");
831 if (token != NULL) extrapH[i] = atoi(token);
832 }
833 } else if (strcmp(token, "L") == 0) {
834 token = strtok(NULL, " :,%=\t");
835 for (int i = 0; i < 13; i++) {
836 token = strtok(NULL, " :,%=\t");
837 if (token != NULL) extrapL[i] = atoi(token);
838 }
839 } else if (strcmp(token, "Thresholds") == 0) {
840 token = strtok(NULL, " :,%=\t");
841 if (token != NULL) thrAlp = atoi(token);
842 token = strtok(NULL, " :,%=\t");
843 if (token != NULL) thrAtt = atoi(token);
844 token = strtok(NULL, " :,%=\t");
845 if (token != NULL) thrDis = atoi(token);
846 } else if (strcmp(token, "Interp") == 0) {
847 for (int i = 0; i < 13; i++) {
848 token = strtok(NULL, " :,%=\t");
849 if (token != NULL) interp[i] = atoi(token);
850 }
851 } else if (strcmp(token, "A") == 0) {
852 // Parameter for energy loss distribution, not used in Garfield++.
853 token = strtok(NULL, " :,%=\t");
854 } else if (strcmp(token, "Z") == 0) {
855 // Parameter for energy loss distribution, not used in Garfield++.
856 token = strtok(NULL, " :,%=\t");
857 } else if (strcmp(token, "EMPROB") == 0) {
858 // Parameter for energy loss distribution, not used in Garfield++.
859 token = strtok(NULL, " :,%=\t");
860 } else if (strcmp(token, "EPAIR") == 0) {
861 // Parameter for energy loss distribution, not used in Garfield++.
862 token = strtok(NULL, " :,%=\t");
863 } else if (strcmp(token, "Ion") == 0) {
864 // Ion diffusion coefficients
865 token = strtok(NULL, " :,%=\t");
866 token = strtok(NULL, " :,%=\t");
867 if (token != NULL) ionDiffL = atof(token);
868 token = strtok(NULL, " :,%=\t");
869 if (token != NULL) ionDiffT = atof(token);
870 } else if (strcmp(token, "CMEAN") == 0) {
871 // Clusters per cm, not used in Garfield..
872 token = strtok(NULL, " :,%=\t");
873 } else if (strcmp(token, "RHO") == 0) {
874 // Parameter for energy loss distribution, not used in Garfield++.
875 token = strtok(NULL, " :,%=\t");
876 } else if (strcmp(token, "PGAS") == 0) {
877 // Pressure [Torr]
878 token = strtok(NULL, " :,%=\t");
879 if (token != NULL) pgas = atof(token);
880 } else if (strcmp(token, "TGAS") == 0) {
881 // Temperature [K]
882 token = strtok(NULL, " :,%=\t");
883 if (token != NULL) tgas = atof(token);
884 done = true;
885 break;
886 } else {
887 done = true;
888 break;
889 }
890 token = strtok(NULL, " :,%=\t");
891 }
892 }
893}
894
895bool MediumGas::GetMixture(const std::vector<double>& mixture,
896 const int version, std::vector<std::string>& gasnames,
897 std::vector<double>& percentages) const {
898
899 gasnames.clear();
900 percentages.clear();
901 const unsigned int nMagboltzGases = mixture.size();
902 for (unsigned int i = 0; i < nMagboltzGases; ++i) {
903 if (mixture[i] < Small) continue;
904 const std::string gasname = GetGasName(i + 1, version);
905 if (gasname.empty()) {
906 std::cerr << m_className << "::GetMixture:\n"
907 << " Unknown gas (gas number " << i + 1 << ").\n";
908 return false;
909 }
910 gasnames.push_back(gasname);
911 percentages.push_back(mixture[i]);
912 }
913 if (gasnames.size() > m_nMaxGases) {
914 std::cerr << m_className << "::GetMixture:\n"
915 << " Gas mixture has " << gasnames.size() << " components.\n"
916 << " Number of gases is limited to " << m_nMaxGases << ".\n";
917 return false;
918 } else if (gasnames.empty()) {
919 std::cerr << m_className << "::GetMixture:\n"
920 << " Gas mixture is not defined (zero components).\n";
921 return false;
922 }
923 double sum = std::accumulate(percentages.begin(), percentages.end(), 0.);
924 if (sum != 100.) {
925 std::cout << m_className << "::GetMixture:\n"
926 << " Renormalizing the percentages.\n";
927 for (auto& percentage : percentages) percentage *= 100. / sum;
928 }
929 return true;
930}
931
932bool MediumGas::MergeGasFile(const std::string& filename,
933 const bool replaceOld) {
934
935 // -----------------------------------------------------------------------
936 // GASMRG - Merges gas data from a file with existing gas tables.
937 // (Last changed on 16/ 2/11.)
938 // -----------------------------------------------------------------------
939
940 constexpr double eps = 1.e-3;
941
942 std::ifstream gasfile;
943 // Open the file.
944 gasfile.open(filename.c_str());
945 // Make sure the file could be opened.
946 if (!gasfile.is_open()) {
947 std::cerr << m_className << "::MergeGasFile:\n"
948 << " Cannot open file " << filename << ".\n";
949 return false;
950 }
951
952 int version = 12;
953 std::bitset<20> gasok;
954 bool new3d = false;
955 constexpr int nMagboltzGases = 60;
956 std::vector<double> mixture(nMagboltzGases, 0.);
957 std::vector<double> efields;
958 std::vector<double> bfields;
959 std::vector<double> angles;
960 std::vector<ExcLevel> excLevels;
961 std::vector<IonLevel> ionLevels;
962 if (!ReadHeader(gasfile, version, gasok, new3d, mixture, efields, bfields,
963 angles, excLevels, ionLevels)) {
964 std::cerr << m_className << "::MergeGasFile: Error reading header.\n";
965 gasfile.close();
966 return false;
967 }
968 // Check the version.
969 if (version != 12) {
970 std::cout << m_className << "::MergeGasFile:\n "
971 << "This dataset cannot be read because of a change in format.\n";
972 gasfile.close();
973 return false;
974 }
975
976 // Check the gas composition.
977 std::vector<std::string> gasnames;
978 std::vector<double> percentages;
979 if (!GetMixture(mixture, version, gasnames, percentages)) {
980 std::cerr << m_className << "::MergeGasFile:\n "
981 << "Cannot determine the gas composition.\n";
982 gasfile.close();
983 return false;
984 }
985 if (m_nComponents != gasnames.size()) {
986 std::cerr << m_className << "::MergeGasFile:\n "
987 << "Composition of the dataset differs from the present one.\n";
988 gasfile.close();
989 return false;
990 }
991
992 for (unsigned int i = 0; i < m_nComponents; ++i) {
993 const auto it = std::find(gasnames.begin(), gasnames.end(), m_gas[i]);
994 if (it == gasnames.end()) {
995 std::cerr << m_className << "::MergeGasFile:\n "
996 << "Composition of the dataset differs from the present one.\n";
997 gasfile.close();
998 return false;
999 }
1000 const double f2 = m_fraction[i];
1001 const double f1 = 0.01 * percentages[it - gasnames.begin()];
1002 if (fabs(f1 - f2) > 1.e-6 * (1. + fabs(f1) + fabs(f2))) {
1003 std::cerr << m_className << "::MergeGasFile:\n "
1004 << "Percentages of " << m_gas[i] << " differ.\n";
1005 gasfile.close();
1006 return false;
1007 }
1008 }
1009
1010 // Check that the excitations and ionisations match.
1011 const unsigned int nexc = excLevels.size();
1012 const unsigned int nion = ionLevels.size();
1013 bool excMatch = (m_excLevels.size() == nexc);
1014 if (excMatch) {
1015 for (unsigned int i = 0; i < nexc; ++i) {
1016 if (m_excLevels[i].label == excLevels[i].label) continue;
1017 excMatch = false;
1018 break;
1019 }
1020 }
1021 bool ionMatch = (m_ionLevels.size() == nion);
1022 if (ionMatch) {
1023 for (unsigned int i = 0; i < nion; ++i) {
1024 if (m_ionLevels[i].label == ionLevels[i].label) continue;
1025 ionMatch = false;
1026 break;
1027 }
1028 }
1029
1030 // Drift velocity along E, Bt and ExB
1031 double ve = 0., vb = 0., vx = 0.;
1032 // Lorentz angle
1033 double lor = 0.;
1034 // Longitudinal and transverse diffusion coefficients
1035 double dl = 0., dt = 0.;
1036 // Townsend and attachment coefficients
1037 double alpha = 0., alpha0 = 0., eta = 0.;
1038 // Ion mobility and dissociation coefficient
1039 double mu = 0., dis = 0.;
1040 // Diffusion tensor.
1041 std::array<double, 6> diff;
1042 // Excitation and ionization rates.
1043 std::vector<double> rexc(nexc, 0.);
1044 std::vector<double> rion(nion, 0.);
1045 // Loop through the gas tables to fast-forward to the footer.
1046 const unsigned int nNewE = efields.size();
1047 const unsigned int nNewB = bfields.size();
1048 const unsigned int nNewA = angles.size();
1049 for (unsigned int i = 0; i < nNewE; i++) {
1050 for (unsigned int j = 0; j < nNewA; j++) {
1051 for (unsigned int k = 0; k < nNewB; k++) {
1052 if (new3d) {
1053 ReadRecord3D(gasfile, ve, vb, vx, dl, dt, alpha, alpha0, eta, mu,
1054 lor, dis, diff, rexc, rion);
1055 } else {
1056 ReadRecord1D(gasfile, ve, vb, vx, dl, dt, alpha, alpha0, eta, mu,
1057 lor, dis, diff, rexc, rion);
1058 }
1059 }
1060 }
1061 }
1062
1063 // Extrapolation methods
1064 std::array<unsigned int, 13> extrapH = {{0}};
1065 std::array<unsigned int, 13> extrapL = {{1}};
1066 // Interpolation methods
1067 std::array<unsigned int, 13> interp = {{2}};
1068 // Thresholds.
1069 unsigned int thrAlp = 0, thrAtt = 0, thrDis = 0;
1070 // Ion diffusion coefficients.
1071 double ionDiffL = 0., ionDiffT = 0.;
1072 // Gas pressure [Torr] and temperature [K].
1073 double pgas = 0., tgas = 0.;
1074 // Read the footer.
1075 gasfile.ignore(std::numeric_limits<std::streamsize>::max(), '\n');
1076 ReadFooter(gasfile, extrapH, extrapL, interp, thrAlp, thrAtt, thrDis,
1077 ionDiffL, ionDiffT, pgas, tgas);
1078
1079 // Check the pressure and temperature.
1080 if (!Similar(pgas, m_pressureTable, eps)) {
1081 std::cerr << m_className << "::MergeGasFile:\n "
1082 << "The gas pressure of the dataset to be read differs\n "
1083 << "from the current reference pressure; stop.\n";
1084 gasfile.close();
1085 return false;
1086 }
1087 if (!Similar(tgas, m_temperatureTable, eps)) {
1088 std::cerr << m_className << "::MergeGasFile:\n "
1089 << "The gas temperature of the dataset to be read differs\n "
1090 << "from the current reference temperature; stop.\n";
1091 gasfile.close();
1092 return false;
1093 }
1094
1095 // Go back to the start and re-read the header to get to the gas tables.
1096 gasfile.clear();
1097 gasfile.seekg(0);
1098 if (!ReadHeader(gasfile, version, gasok, new3d, mixture, efields, bfields,
1099 angles, excLevels, ionLevels)) {
1100 std::cerr << m_className << "::MergeGasFile: Error re-reading header.\n";
1101 gasfile.close();
1102 return false;
1103 }
1104
1105 // Multiply the E/p values by the pressure.
1106 for (auto& field : efields) field *= pgas;
1107
1108 if (m_debug) {
1109 std::cout << m_className << "::MergeGasFile:\n "
1110 << "Dataset to be merged has the following dimensions:\n "
1111 << "3D = " << new3d << " nE = " << nNewE << ", nB = " << nNewB
1112 << ", nA = " << nNewA << ", nExc = "
1113 << excLevels.size() << ", nIon = " << ionLevels.size() << "\n";
1114 }
1115
1116 unsigned int nE = m_eFields.size();
1117 unsigned int nB = m_bFields.size();
1118 unsigned int nA = m_bAngles.size();
1119 // Determine which mode we have to use for the E field.
1120 const int iemode = Equal(efields, m_eFields, eps);
1121 // Determine which mode we have to use for the angles.
1122 const int iamode = Equal(angles, m_bAngles, eps);
1123 // Determine which mode we have to use for the B field.
1124 const int ibmode = Equal(bfields, m_bFields, eps);
1125 if (m_debug) {
1126 std::cout << m_className << "::MergeGasFile:\n";
1127 if (iemode == 0) std::cout << " The E vectors differ.\n";
1128 else std::cout << " The E vectors are identical.\n";
1129 if (iamode == 0) std::cout << " The angle vectors differ.\n";
1130 else std::cout << " The angle vectors are identical.\n";
1131 if (ibmode == 0) std::cout << " The B vectors differ.\n";
1132 else std::cout << " The B vectors are identical.\n";
1133 }
1134 // Ensure there is a common mode.
1135 if (iemode + iamode + ibmode < 2) {
1136 std::cerr << m_className << "::MergeGasFile:\n Existing data and data "
1137 << "in the file don't have two common axes; not merged.\n";
1138 gasfile.close();
1139 return false;
1140 }
1141 // Decide whether we have to produce a 3D table or a 1D table.
1142 const bool old3d = m_tab2d;
1143 if ((new3d || ibmode * iamode == 0) && !m_tab2d) {
1144 if (m_debug) std::cout << " Expanding existing table to 3D mode.\n";
1145 m_tab2d = true;
1146 }
1147
1148 // Determine which data are currently present.
1149 std::bitset<20> existing;
1150 GetGasBits(existing);
1151 // If the grids don't match, warn for data being lost in the merge.
1152 if (iemode * ibmode * iamode == 0) {
1153 // Check for data currently present which are absent in the new data.
1154 if (existing[0] && !gasok[0]) {
1155 existing.reset(0);
1156 m_eVelE.clear();
1157 PrintAbsentInNew("drift velocity");
1158 }
1159 if (existing[1] && !gasok[1]) {
1160 existing.reset(1);
1161 m_iMob.clear();
1162 PrintAbsentInNew("ion mobility");
1163 }
1164 if (existing[2] && !gasok[2]) {
1165 existing.reset(2);
1166 m_eDifL.clear();
1167 PrintAbsentInNew("longitudinal diffusion");
1168 }
1169 if (existing[3] && !gasok[3]) {
1170 existing.reset(3);
1171 m_eAlp.clear();
1172 PrintAbsentInNew("Townsend coefficient");
1173 }
1174 if (existing[5] && !gasok[5]) {
1175 existing.reset(5);
1176 m_eAtt.clear();
1177 PrintAbsentInNew("attachment coefficient");
1178 }
1179 if (existing[6] && !gasok[6]) {
1180 existing.reset(6);
1181 m_eLor.clear();
1182 PrintAbsentInNew("Lorentz angle");
1183 }
1184 if (existing[7] && !gasok[7]) {
1185 existing.reset(7);
1186 m_eDifT.clear();
1187 PrintAbsentInNew("transverse diffusion");
1188 }
1189 if (existing[8] && !gasok[8]) {
1190 existing.reset(8);
1191 m_eVelB.clear();
1192 PrintAbsentInNew("velocity along Bt");
1193 }
1194 if (existing[9] && !gasok[9]) {
1195 existing.reset(9);
1196 m_eVelX.clear();
1197 PrintAbsentInNew("velocity along ExB");
1198 }
1199 if (existing[10] && !gasok[10]) {
1200 existing.reset(10);
1201 m_eDifM.clear();
1202 PrintAbsentInNew("diffusion tensor");
1203 }
1204 if (existing[11] && !gasok[11]) {
1205 existing.reset(11);
1206 m_iDis.clear();
1207 PrintAbsentInNew("ion dissociation data");
1208 }
1209 if (existing[14] && !gasok[14]) {
1210 existing.reset(14);
1211 m_excLevels.clear();
1212 m_excRates.clear();
1213 PrintAbsentInNew("excitation data");
1214 }
1215 if (existing[15] && !gasok[15]) {
1216 existing.reset(15);
1217 m_ionLevels.clear();
1218 m_ionRates.clear();
1219 PrintAbsentInNew("ionisation data");
1220 }
1221 // And for data present in the file but not currently present.
1222 if (!existing[0] && gasok[0]) {
1223 gasok.reset(0);
1224 PrintAbsentInExisting("drift velocity");
1225 }
1226 if (!existing[1] && gasok[1]) {
1227 gasok.reset(1);
1228 PrintAbsentInExisting("ion mobility");
1229 }
1230 if (!existing[2] && gasok[2]) {
1231 gasok.reset(2);
1232 PrintAbsentInExisting("longitudinal diffusion");
1233 }
1234 if (!existing[3] && gasok[3]) {
1235 gasok.reset(3);
1236 PrintAbsentInExisting("Townsend coefficient");
1237 }
1238 if (!existing[5] && gasok[5]) {
1239 gasok.reset(5);
1240 PrintAbsentInExisting("attachment coefficient");
1241 }
1242 if (!existing[6] && gasok[6]) {
1243 if (old3d) {
1244 gasok.reset(6);
1245 PrintAbsentInExisting("Lorentz angle");
1246 } else {
1247 // Existing dataset is for B = 0 (no Lorentz angle).
1248 Init(nE, nB, nA, m_eLor, -30.);
1249 existing.set(6);
1250 }
1251 }
1252 if (!existing[7] && gasok[7]) {
1253 gasok.reset(7);
1254 PrintAbsentInExisting("transverse diffusion");
1255 }
1256 if (!existing[8] && gasok[8]) {
1257 if (old3d) {
1258 gasok.reset(8);
1259 PrintAbsentInExisting("velocity along Bt");
1260 } else {
1261 // Existing dataset is for B = 0 (no velocity component || Bt).
1262 Init(nE, nB, nA, m_eVelB, 0.);
1263 existing.set(8);
1264 }
1265 }
1266 if (!existing[9] && gasok[9]) {
1267 if (old3d) {
1268 gasok.reset(9);
1269 PrintAbsentInExisting("velocity along ExB");
1270 } else {
1271 // Existing dataset is for B = 0 (no velocity component || ExB).
1272 Init(nE, nB, nA, m_eVelX, 0.);
1273 existing.set(9);
1274 }
1275 }
1276 if (!existing[10] && gasok[10]) {
1277 gasok.reset(10);
1278 PrintAbsentInExisting("diffusion tensor");
1279 }
1280 if (!existing[11] && gasok[11]) {
1281 gasok.reset(11);
1282 PrintAbsentInExisting("ion dissociation data");
1283 }
1284 if (!existing[14] && gasok[14]) {
1285 gasok.reset(14);
1286 PrintAbsentInExisting("excitation data");
1287 }
1288 if (!existing[15] && gasok[15]) {
1289 gasok.reset(15);
1290 PrintAbsentInExisting("ionisation data");
1291 }
1292 if (existing[14] && gasok[14] && !excMatch) {
1293 std::cerr << " Excitation levels of the two datasets don't match.\n"
1294 << " Deleting excitation data.\n";
1295 m_excLevels.clear();
1296 m_excRates.clear();
1297 existing.reset(14);
1298 gasok.reset(14);
1299 }
1300 if (existing[15] && gasok[15] && !ionMatch) {
1301 std::cerr << " Ionisation levels of the two datasets don't match.\n"
1302 << " Deleting ionisation data.\n";
1303 m_ionLevels.clear();
1304 m_ionRates.clear();
1305 existing.reset(15);
1306 gasok.reset(15);
1307 }
1308 } else {
1309 // If the grids are identical, initialise the tables that are only present
1310 // in the new dataset but not in existing one.
1311 if (gasok[0] && !existing[0]) Init(nE, nB, nA, m_eVelE, 0.);
1312 if (gasok[1] && !existing[1]) Init(nE, nB, nA, m_iMob, 0.);
1313 if (gasok[2] && !existing[2]) Init(nE, nB, nA, m_eDifL, 0.);
1314 if (gasok[3] && !existing[3]) {
1315 Init(nE, nB, nA, m_eAlp, -30.);
1316 Init(nE, nB, nA, m_eAlp0, -30.);
1317 }
1318 if (gasok[5] && !existing[5]) Init(nE, nB, nA, m_eAtt, -30.);
1319 if (gasok[6] && !existing[6]) Init(nE, nB, nA, m_eLor, -30.);
1320 if (gasok[7] && !existing[7]) Init(nE, nB, nA, m_eDifT, 0.);
1321 if (gasok[8] && !existing[8]) Init(nE, nB, nA, m_eVelB, 0.);
1322 if (gasok[9] && !existing[9]) Init(nE, nB, nA, m_eVelX, 0.);
1323 if (gasok[10] && !existing[10]) Init(nE, nB, nA, 6, m_eDifM, 0.);
1324 if (gasok[11] && !existing[11]) Init(nE, nB, nA, m_iDis, -30.);
1325 if (gasok[14] && (!existing[14] || replaceOld)) {
1326 Init(nE, nB, nA, nexc, m_excRates, 0.);
1327 }
1328 if (gasok[15] && (!existing[15] || replaceOld)) {
1329 Init(nE, nB, nA, nion, m_ionRates, 0.);
1330 }
1331 }
1332
1333 // Initialise the "new" flags.
1334 std::vector<bool> newE(nE, false);
1335 std::vector<bool> newA(nA, false);
1336 std::vector<bool> newB(nB, false);
1337 // Extend the existing tables.
1338 std::cout << m_className << "::MergeGasFile: Extending the tables.\n";
1339 // Insert room in the tables for new columns in E.
1340 if (iemode == 0) {
1341 // Loop over the new values.
1342 for (const auto efield : efields) {
1343 // Loop over the old values.
1344 bool found = false;
1345 for (unsigned int j = 0; j < nE; ++j) {
1346 // If it overlaps with existing E, either keep old or new data.
1347 if (Similar(efield, m_eFields[j], eps)) {
1348 if (replaceOld) {
1349 std::cout << " Replacing existing data for E = "
1350 << m_eFields[j] << " V/cm by data from file.\n";
1351 m_eFields[j] = efield;
1352 newE[j] = true;
1353 ZeroRowE(j, nB, nA);
1354 } else {
1355 std::cout << " Keeping existing data for E = " << m_eFields[j]
1356 << " V/cm, not using data from the file.\n";
1357 }
1358 found = true;
1359 break;
1360 } else if (efield < m_eFields[j]) {
1361 // Otherwise shift all data at higher E values.
1362 if (m_debug) {
1363 std::cout << " Inserting E = " << efield
1364 << " V/cm at slot " << j << ".\n";
1365 }
1366 InsertE(j, nE, nB, nA);
1367 m_eFields.insert(m_eFields.begin() + j, efield);
1368 newE.insert(newE.begin() + j, true);
1369 ZeroRowE(j, nB, nA);
1370 ++nE;
1371 found = true;
1372 break;
1373 }
1374 }
1375 if (found) continue;
1376 // If there is no higher E, then add the line at the end.
1377 if (m_debug) {
1378 std::cout << " Adding E = " << efield << " V/cm at the end.\n";
1379 }
1380 InsertE(nE, nE, nB, nA);
1381 m_eFields.push_back(efield);
1382 newE.push_back(true);
1383 ZeroRowE(nE, nB, nA);
1384 ++nE;
1385 }
1386 }
1387 // Insert room in the tables for new columns in B.
1388 if (ibmode == 0) {
1389 // Loop over the new values.
1390 for (const auto bfield : bfields) {
1391 // Loop over the old values.
1392 bool found = false;
1393 for (unsigned int j = 0; j < nB; ++j) {
1394 // If it overlaps with existing B, either keep old or new data.
1395 if (Similar(bfield, m_bFields[j], eps)) {
1396 if (replaceOld) {
1397 std::cout << " Replacing old data for B = " << m_bFields[j]
1398 << " T by data from file.\n";
1399 m_bFields[j] = bfield;
1400 newB[j] = true;
1401 ZeroRowB(j, nE, nA);
1402 } else {
1403 std::cout << " Keeping old data for B = " << m_bFields[j]
1404 << " T, not using data from file.\n";
1405 }
1406 found = true;
1407 break;
1408 } else if (bfield < m_bFields[j]) {
1409 // Otherwise shift all data at higher B values.
1410 if (m_debug) {
1411 std::cout << " Inserting B = " << bfield << " T at slot "
1412 << j << ".\n";
1413 }
1414 InsertB(j, nE, nB, nA);
1415 m_bFields.insert(m_bFields.begin() + j, bfield);
1416 newB.insert(newB.begin() + j, true);
1417 ZeroRowB(j, nE, nA);
1418 ++nB;
1419 found = true;
1420 break;
1421 }
1422 }
1423 if (found) continue;
1424 // If there is no higher B, then add the line at the end.
1425 if (m_debug) {
1426 std::cout << " Adding B = " << bfield << " T at the end.\n";
1427 }
1428 InsertB(nB, nE, nB, nA);
1429 m_bFields.push_back(bfield);
1430 newB.push_back(true);
1431 ZeroRowB(nB, nE, nA);
1432 ++nB;
1433 }
1434 }
1435 // Insert room in the tables for new columns in angle.
1436 if (iamode == 0) {
1437 // Loop over the new values.
1438 for (const auto angle : angles) {
1439 // Loop over the old values.
1440 bool found = false;
1441 for (unsigned int j = 0; j < nA; ++j) {
1442 // If it overlaps with an existing angle, either keep old or new data.
1443 if (Similar(angle, m_bAngles[j], eps)) {
1444 if (replaceOld) {
1445 std::cout << " Replacing old data for angle(E,B) = "
1446 << m_bAngles[j] * RadToDegree
1447 << " degrees by data from the file.\n";
1448 m_bAngles[j] = angle;
1449 newA[j] = true;
1450 ZeroRowA(j, nE, nB);
1451 } else {
1452 std::cout << " Keeping old data for angle(E,B) = "
1453 << m_bAngles[j] * RadToDegree
1454 << " degrees, not using data from file.\n";
1455 }
1456 found = true;
1457 break;
1458 } else if (angle < m_bAngles[j]) {
1459 // Otherwise shift all data at higher angles.
1460 if (m_debug) {
1461 std::cout << " Inserting angle = " << angle * RadToDegree
1462 << " degrees at slot " << j << ".\n";
1463 }
1464 InsertA(j, nE, nB, nA);
1465 m_bAngles.insert(m_bAngles.begin() + j, angle);
1466 newA.insert(newA.begin() + j, true);
1467 ZeroRowA(j, nE, nB);
1468 ++nA;
1469 found = true;
1470 break;
1471 }
1472 }
1473 if (found) continue;
1474 // If there is no higher angle, then add the line at the end.
1475 if (m_debug) {
1476 std::cout << " Adding angle = " << angle * RadToDegree
1477 << " degrees at the end.\n";
1478 }
1479 InsertA(nA, nE, nB, nA);
1480 m_bAngles.push_back(angle);
1481 newA.push_back(true);
1482 ZeroRowA(nA, nE, nB);
1483 ++nA;
1484 }
1485 }
1486
1487 const double sqrp = sqrt(pgas);
1488 const double logp = log(pgas);
1489
1490 // Read the gas table.
1491 for (const auto efield : efields) {
1492 // Locate the index at which these values are to be stored.
1493 const int inde = FindIndex(efield, m_eFields, eps);
1494 for (const auto angle : angles) {
1495 const int inda = FindIndex(angle, m_bAngles, eps);
1496 for (const auto bfield : bfields) {
1497 // Read the record.
1498 if (new3d) {
1499 ReadRecord3D(gasfile, ve, vb, vx, dl, dt, alpha, alpha0, eta, mu,
1500 lor, dis, diff, rexc, rion);
1501 } else {
1502 ReadRecord1D(gasfile, ve, vb, vx, dl, dt, alpha, alpha0, eta, mu,
1503 lor, dis, diff, rexc, rion);
1504 }
1505 const int indb = FindIndex(bfield, m_bFields, eps);
1506 if (inde < 0 || inda < 0 || indb < 0) {
1507 std::cerr << m_className << "::MergeGasFile:\n Unable to locate"
1508 << " the (E,angle,B) insertion point; no gas data read.\n";
1509 std::cout << "BFIELD = " << bfield << ", IB = " << indb << "\n";
1510 ResetTables();
1511 gasfile.close();
1512 return false;
1513 }
1514 const bool update = newE[inde] || newA[inda] || newB[indb] || replaceOld;
1515 // Store the data.
1516 if (gasok[0] && (update || !existing[0])) {
1517 m_eVelE[inda][indb][inde] = ve;
1518 }
1519 if (gasok[1] && (update || !existing[1])) {
1520 m_iMob[inda][indb][inde] = mu;
1521 }
1522 if (gasok[2] && (update || !existing[2])) {
1523 m_eDifL[inda][indb][inde] = dl / sqrp;
1524 }
1525 if (gasok[3] && (update || !existing[3])) {
1526 m_eAlp[inda][indb][inde] = alpha + logp;
1527 m_eAlp0[inda][indb][inde] = alpha0 + logp;
1528 }
1529 if (gasok[5] && (update || !existing[5])) {
1530 m_eAtt[inda][indb][inde] = eta + logp;
1531 }
1532 if (gasok[6] && (update || !existing[6])) {
1533 m_eLor[inda][indb][inde] = lor;
1534 }
1535 if (gasok[7] && (update || !existing[7])) {
1536 m_eDifT[inda][indb][inde] = dt / sqrp;
1537 }
1538 if (gasok[8] && (update || !existing[8])) {
1539 m_eVelB[inda][indb][inde] = vb;
1540 }
1541 if (gasok[9] && (update || !existing[9])) {
1542 m_eVelX[inda][indb][inde] = vx;
1543 }
1544 if (gasok[10] && (update || !existing[10])) {
1545 for (unsigned int l = 0; l < 6; ++l) {
1546 m_eDifM[l][inda][indb][inde] = diff[l] / pgas;
1547 }
1548 }
1549 if (gasok[11] && (update || !existing[11])) {
1550 m_iDis[inda][indb][inde] = dis + logp;
1551 }
1552 if (gasok[14] && (update || !existing[14])) {
1553 for (unsigned int l = 0; l < nexc; ++l) {
1554 m_excRates[l][inda][indb][inde] = rexc[l];
1555 }
1556 }
1557 if (gasok[15] && (update || !existing[15])) {
1558 for (unsigned int l = 0; l < nion; ++l) {
1559 m_ionRates[l][inda][indb][inde] = rion[l];
1560 }
1561 }
1562 }
1563 }
1564 }
1565 // if (iemode + iamode + ibmode == 3) { ... }
1566 if (replaceOld) {
1567 if (m_debug) {
1568 std::cout << m_className << "::MergeGasFile: "
1569 << "Replacing extrapolation and interpolation data.\n";
1570 }
1571 if (gasok[0]) m_extrVel = {extrapL[0], extrapH[0]};
1572 if (gasok[1]) m_extrMob = {extrapL[6], extrapH[6]};
1573 if (gasok[2]) m_extrDif = {extrapL[3], extrapH[3]};
1574 if (gasok[3]) m_extrAlp = {extrapL[4], extrapH[4]};
1575 if (gasok[5]) m_extrAtt = {extrapL[5], extrapH[5]};
1576 if (gasok[6]) m_extrLor = {extrapL[7], extrapH[7]};
1577 if (gasok[11]) m_extrDis = {extrapL[9], extrapH[9]};
1578 if (gasok[14]) m_extrExc = {extrapL[11], extrapH[11]};
1579 if (gasok[15]) m_extrIon = {extrapL[12], extrapH[12]};
1580
1581 if (gasok[0]) m_intpVel = interp[0];
1582 if (gasok[1]) m_intpMob = interp[6];
1583 if (gasok[2]) m_intpDif = interp[3];
1584 if (gasok[3]) m_intpAlp = interp[4];
1585 if (gasok[5]) m_intpAtt = interp[5];
1586 if (gasok[6]) m_intpLor = interp[7];
1587 if (gasok[11]) m_intpDis = interp[9];
1588 if (gasok[14]) m_intpExc = interp[11];
1589 if (gasok[15]) m_intpIon = interp[12];
1590
1591 // Ion diffusion.
1592 if (m_debug && (ionDiffL > 0. || ionDiffT > 0.)) {
1593 std::cout << m_className << "::MergeGasFile: Replacing ion diffusion.\n";
1594 }
1595 if (ionDiffL > 0.) Init(nE, nB, nA, m_iDifL, ionDiffL);
1596 if (ionDiffT > 0.) Init(nE, nB, nA, m_iDifT, ionDiffT);
1597 }
1598 // Update the Townsend and attachment threshold indices.
1601 return true;
1602}
1603
1604void MediumGas::InsertE(const int ie, const int ne, const int nb,
1605 const int na) {
1606 for (int k = 0; k < na; ++k) {
1607 for (int j = 0; j < nb; ++j) {
1608 if (!m_eVelE.empty()) m_eVelE[k][j].resize(ne + 1, 0.);
1609 if (!m_eVelB.empty()) m_eVelB[k][j].resize(ne + 1, 0.);
1610 if (!m_eVelX.empty()) m_eVelX[k][j].resize(ne + 1, 0.);
1611 if (!m_eDifL.empty()) m_eDifL[k][j].resize(ne + 1, 0.);
1612 if (!m_eDifT.empty()) m_eDifT[k][j].resize(ne + 1, 0.);
1613 if (!m_eAlp.empty()) m_eAlp[k][j].resize(ne + 1, 0.);
1614 if (!m_eAlp0.empty()) m_eAlp0[k][j].resize(ne + 1, 0.);
1615 if (!m_eAtt.empty()) m_eAtt[k][j].resize(ne + 1, 0.);
1616 if (!m_eLor.empty()) m_eLor[k][j].resize(ne + 1, 0.);
1617 if (!m_iMob.empty()) m_iMob[k][j].resize(ne + 1, 0.);
1618 if (!m_iDis.empty()) m_iDis[k][j].resize(ne + 1, 0.);
1619 if (!m_iDifL.empty()) m_iDifL[k][j].resize(ne + 1, 0.);
1620 if (!m_iDifT.empty()) m_iDifT[k][j].resize(ne + 1, 0.);
1621 for (auto& dif : m_eDifM) dif[k][j].resize(ne + 1, 0.);
1622 for (auto& exc : m_excRates) exc[k][j].resize(ne + 1, 0.);
1623 for (auto& ion : m_ionRates) ion[k][j].resize(ne + 1, 0.);
1624 for (int i = ne; i > ie; --i) {
1625 if (!m_eVelE.empty()) m_eVelE[k][j][i] = m_eVelE[k][j][i - 1];
1626 if (!m_eVelB.empty()) m_eVelB[k][j][i] = m_eVelB[k][j][i - 1];
1627 if (!m_eVelX.empty()) m_eVelX[k][j][i] = m_eVelX[k][j][i - 1];
1628 if (!m_eDifL.empty()) m_eDifL[k][j][i] = m_eDifL[k][j][i - 1];
1629 if (!m_eDifT.empty()) m_eDifT[k][j][i] = m_eDifT[k][j][i - 1];
1630 if (!m_eAlp.empty()) m_eAlp[k][j][i] = m_eAlp[k][j][i - 1];
1631 if (!m_eAlp0.empty()) m_eAlp0[k][j][i] = m_eAlp0[k][j][i - 1];
1632 if (!m_eAtt.empty()) m_eAtt[k][j][i] = m_eAtt[k][j][i - 1];
1633 if (!m_eLor.empty()) m_eLor[k][j][i] = m_eLor[k][j][i - 1];
1634 if (!m_iMob.empty()) m_iMob[k][j][i] = m_iMob[k][j][i - 1];
1635 if (!m_iDis.empty()) m_iDis[k][j][i] = m_iDis[k][j][i - 1];
1636 if (!m_iDifL.empty()) m_iDifL[k][j][i] = m_iDifL[k][j][i - 1];
1637 if (!m_iDifT.empty()) m_iDifT[k][j][i] = m_iDifT[k][j][i - 1];
1638 for (auto& dif : m_eDifM) dif[k][j][i] = dif[k][j][i - 1];
1639 for (auto& exc : m_excRates) exc[k][j][i] = exc[k][j][i - 1];
1640 for (auto& ion : m_ionRates) ion[k][j][i] = ion[k][j][i - 1];
1641 }
1642 }
1643 }
1644}
1645
1646void MediumGas::InsertB(const int ib, const int ne, const int nb,
1647 const int na) {
1648 for (int k = 0; k < na; ++k) {
1649 if (!m_eVelE.empty()) m_eVelE[k].resize(nb + 1, std::vector<double>(ne, 0.));
1650 if (!m_eVelB.empty()) m_eVelB[k].resize(nb + 1, std::vector<double>(ne, 0.));
1651 if (!m_eVelX.empty()) m_eVelX[k].resize(nb + 1, std::vector<double>(ne, 0.));
1652 if (!m_eDifL.empty()) m_eDifL[k].resize(nb + 1, std::vector<double>(ne, 0.));
1653 if (!m_eDifT.empty()) m_eDifT[k].resize(nb + 1, std::vector<double>(ne, 0.));
1654 if (!m_eAlp.empty()) m_eAlp[k].resize(nb + 1, std::vector<double>(ne, 0.));
1655 if (!m_eAlp0.empty()) m_eAlp0[k].resize(nb + 1, std::vector<double>(ne, 0.));
1656 if (!m_eAtt.empty()) m_eAtt[k].resize(nb + 1, std::vector<double>(ne, 0.));
1657 if (!m_eLor.empty()) m_eLor[k].resize(nb + 1, std::vector<double>(ne, 0.));
1658 if (!m_iMob.empty()) m_iMob[k].resize(nb + 1, std::vector<double>(ne, 0.));
1659 if (!m_iDis.empty()) m_iDis[k].resize(nb + 1, std::vector<double>(ne, 0.));
1660 if (!m_iDifL.empty()) m_iDifL[k].resize(nb + 1, std::vector<double>(ne, 0.));
1661 if (!m_iDifT.empty()) m_iDifT[k].resize(nb + 1, std::vector<double>(ne, 0.));
1662 for (auto& dif : m_eDifM) {
1663 dif[k].resize(nb + 1, std::vector<double>(ne, 0.));
1664 }
1665 for (auto& exc : m_excRates) {
1666 exc[k].resize(nb + 1, std::vector<double>(ne, 0.));
1667 }
1668 for (auto& ion : m_ionRates) {
1669 ion[k].resize(nb + 1, std::vector<double>(ne, 0.));
1670 }
1671 for (int i = 0; i < ne; ++i) {
1672 for (int j = nb; j > ib; j--) {
1673 if (!m_eVelE.empty()) m_eVelE[k][j][i] = m_eVelE[k][j - 1][i];
1674 if (!m_eVelB.empty()) m_eVelB[k][j][i] = m_eVelB[k][j - 1][i];
1675 if (!m_eVelX.empty()) m_eVelX[k][j][i] = m_eVelX[k][j - 1][i];
1676 if (!m_eDifL.empty()) m_eDifL[k][j][i] = m_eDifL[k][j - 1][i];
1677 if (!m_eDifT.empty()) m_eDifT[k][j][i] = m_eDifT[k][j - 1][i];
1678 if (!m_eAlp.empty()) m_eAlp[k][j][i] = m_eAlp[k][j - 1][i];
1679 if (!m_eAlp0.empty()) m_eAlp0[k][j][i] = m_eAlp0[k][j - 1][i];
1680 if (!m_eAtt.empty()) m_eAtt[k][j][i] = m_eAtt[k][j - 1][i];
1681 if (!m_eLor.empty()) m_eLor[k][j][i] = m_eLor[k][j - 1][i];
1682 if (!m_iMob.empty()) m_iMob[k][j][i] = m_iMob[k][j - 1][i];
1683 if (!m_iDis.empty()) m_iDis[k][j][i] = m_iDis[k][j - 1][i];
1684 if (!m_iDifL.empty()) m_iDifL[k][j][i] = m_iDifL[k][j - 1][i];
1685 if (!m_iDifT.empty()) m_iDifT[k][j][i] = m_iDifT[k][j - 1][i];
1686 for (auto& dif : m_eDifM) dif[k][j][i] = dif[k][j - 1][i];
1687 for (auto& exc : m_excRates) exc[k][j][i] = exc[k][j - 1][i];
1688 for (auto& ion : m_ionRates) ion[k][j][i] = ion[k][j - 1][i];
1689 }
1690 }
1691 }
1692}
1693
1694void MediumGas::InsertA(const int ia, const int ne, const int nb,
1695 const int na) {
1696 ResizeA(m_eVelE, ne, nb, na + 1);
1697 ResizeA(m_eVelB, ne, nb, na + 1);
1698 ResizeA(m_eVelX, ne, nb, na + 1);
1699 ResizeA(m_eDifL, ne, nb, na + 1);
1700 ResizeA(m_eDifT, ne, nb, na + 1);
1701 ResizeA(m_eAlp, ne, nb, na + 1);
1702 ResizeA(m_eAlp0, ne, nb, na + 1);
1703 ResizeA(m_eAtt, ne, nb, na + 1);
1704 ResizeA(m_eLor, ne, nb, na + 1);
1705 ResizeA(m_iMob, ne, nb, na + 1);
1706 ResizeA(m_iDis, ne, nb, na + 1);
1707 ResizeA(m_iDifL, ne, nb, na + 1);
1708 ResizeA(m_iDifT, ne, nb, na + 1);
1709 for (auto& dif : m_eDifM) ResizeA(dif, ne, nb, na + 1);
1710 for (auto& exc : m_excRates) ResizeA(exc, ne, nb, na + 1);
1711 for (auto& ion : m_ionRates) ResizeA(ion, ne, nb, na + 1);
1712 for (int j = 0; j < nb; ++j) {
1713 for (int i = 0; i < ne; ++i) {
1714 for (int k = na; k > ia; k--) {
1715 if (!m_eVelE.empty()) m_eVelE[k][j][i] = m_eVelE[k - 1][j][i];
1716 if (!m_eVelB.empty()) m_eVelB[k][j][i] = m_eVelB[k - 1][j][i];
1717 if (!m_eVelX.empty()) m_eVelX[k][j][i] = m_eVelX[k - 1][j][i];
1718 if (!m_eDifL.empty()) m_eDifL[k][j][i] = m_eDifL[k - 1][j][i];
1719 if (!m_eDifT.empty()) m_eDifT[k][j][i] = m_eDifT[k - 1][j][i];
1720 if (!m_eAlp.empty()) m_eAlp[k][j][i] = m_eAlp[k - 1][j][i];
1721 if (!m_eAlp0.empty()) m_eAlp0[k][j][i] = m_eAlp0[k - 1][j][i];
1722 if (!m_eAtt.empty()) m_eAtt[k][j][i] = m_eAtt[k - 1][j][i];
1723 if (!m_eLor.empty()) m_eLor[k][j][i] = m_eLor[k - 1][j][i];
1724 if (!m_iMob.empty()) m_iMob[k][j][i] = m_iMob[k - 1][j][i];
1725 if (!m_iDis.empty()) m_iDis[k][j][i] = m_iDis[k - 1][j][i];
1726 if (!m_iDifL.empty()) m_iDifL[k][j][i] = m_iDifL[k - 1][j][i];
1727 if (!m_iDifT.empty()) m_iDifT[k][j][i] = m_iDifT[k - 1][j][i];
1728 for (auto& dif : m_eDifM) dif[k][j][i] = dif[k - 1][j][i];
1729 for (auto& exc : m_excRates) exc[k][j][i] = exc[k - 1][j][i];
1730 for (auto& ion : m_ionRates) ion[k][j][i] = ion[k - 1][j][i];
1731 }
1732 }
1733 }
1734}
1735
1736void MediumGas::ZeroRowE(const int ie, const int nb, const int na) {
1737 for (int k = 0; k < na; ++k) {
1738 for (int j = 0; j < nb; ++j) {
1739 if (!m_eVelE.empty()) m_eVelE[k][j][ie] = 0.;
1740 }
1741 }
1742}
1743
1744void MediumGas::ZeroRowB(const int ib, const int ne, const int na) {
1745 for (int k = 0; k < na; ++k) {
1746 for (int i = 0; i < ne; ++i) {
1747 if (!m_eVelE.empty()) m_eVelE[k][ib][i] = 0.;
1748 }
1749 }
1750}
1751
1752void MediumGas::ZeroRowA(const int ia, const int ne, const int nb) {
1753 for (int j = 0; j < nb; ++j) {
1754 for (int i = 0; i < ne; ++i) {
1755 if (!m_eVelE.empty()) m_eVelE[ia][j][i] = 0.;
1756 }
1757 }
1758}
1759
1760bool MediumGas::WriteGasFile(const std::string& filename) {
1761
1762 // -----------------------------------------------------------------------
1763 // GASWRT
1764 // -----------------------------------------------------------------------
1765
1766 // Set the gas mixture.
1767 constexpr int nMagboltzGases = 60;
1768 std::vector<double> mixture(nMagboltzGases, 0.);
1769 for (unsigned int i = 0; i < m_nComponents; ++i) {
1770 const int ng = GetGasNumberGasFile(m_gas[i]);
1771 if (ng <= 0) {
1772 std::cerr << m_className << "::WriteGasFile:\n"
1773 << " Error retrieving gas number for " << m_gas[i] << ".\n";
1774 continue;
1775 }
1776 mixture[ng - 1] = m_fraction[i] * 100.;
1777 }
1778
1779 if (m_debug) {
1780 std::cout << m_className << "::WriteGasFile:\n"
1781 << " Writing gas tables to file " << filename << "\n";
1782 }
1783
1784 std::ofstream outfile;
1785 outfile.open(filename.c_str(), std::ios::out);
1786 if (!outfile.is_open()) {
1787 std::cerr << m_className << "::WriteGasFile:\n"
1788 << " Cannot open file " << filename << ".\n";
1789 outfile.close();
1790 return false;
1791 }
1792
1793 // Assemble the GASOK bits.
1794 std::bitset<20> gasok;
1795 GetGasBits(gasok);
1796 std::string okstr(20, 'F');
1797 for (unsigned int i = 0; i < 20; ++i) {
1798 if (gasok[i]) okstr[i] = 'T';
1799 }
1800 if (m_debug) std::cout << " GASOK bits: " << okstr << "\n";
1801 // Get the current time.
1802 time_t rawtime = time(0);
1803 tm timeinfo = *localtime(&rawtime);
1804 char datebuf[80] = {0};
1805 char timebuf[80] = {0};
1806 // Format date and time.
1807 strftime(datebuf, sizeof(datebuf) - 1, "%d/%m/%y", &timeinfo);
1808 strftime(timebuf, sizeof(timebuf) - 1, "%H.%M.%S", &timeinfo);
1809 // Set the member name.
1810 std::string member = "< none >";
1811 // Write the header.
1812 outfile << "*----.----1----.----2----.----3----.----4----.----"
1813 << "5----.----6----.----7----.----8----.----9----.---"
1814 << "10----.---11----.---12----.---13--\n";
1815 outfile << "% Created " << datebuf << " at " << timebuf << " ";
1816 outfile << member << " GAS ";
1817 // Add remark.
1818 std::string buffer;
1819 outfile << "\"none" << std::string(25, ' ') << "\"\n";
1820 const int version = 12;
1821 outfile << " Version : " << version << "\n";
1822 outfile << " GASOK bits: " << okstr << "\n";
1823 std::stringstream ids;
1824 ids.str("");
1825 for (unsigned int i = 0; i < m_nComponents; ++i) {
1826 ids << m_gas[i] << " " << 100. * m_fraction[i] << "%, ";
1827 }
1828 ids << "T=" << m_temperatureTable << " K, "
1829 << "p=" << m_pressureTable / AtmosphericPressure << " atm";
1830 outfile << " Identifier: " << std::setw(80) << std::left << ids.str() << "\n";
1831 outfile << " Clusters : " << std::string(80, ' ') << "\n";
1832 outfile << " Dimension : ";
1833 if (m_tab2d) {
1834 outfile << "T ";
1835 } else {
1836 outfile << "F ";
1837 }
1838
1839 const unsigned int nE = m_eFields.size();
1840 const unsigned int nB = m_bFields.size();
1841 const unsigned int nA = m_bAngles.size();
1842 if (m_debug) {
1843 std::cout << m_className << "::WriteGasFile:\n "
1844 << "Dataset has the following dimensions:\n "
1845 << "3D = " << m_tab2d << " nE = " << nE << ", nB = " << nB
1846 << ", nA = " << nA << ", nExc = "
1847 << m_excLevels.size() << ", nIon = " << m_ionLevels.size() << "\n";
1848 }
1849 outfile << FmtInt(nE, 9) << " " << FmtInt(nA, 9) << " "
1850 << FmtInt(nB, 9) << " " << FmtInt(m_excLevels.size(), 9) << " "
1851 << FmtInt(m_ionLevels.size(), 9) << "\n";
1852 // Store reduced electric fields (E/p).
1853 outfile << " E fields \n";
1854 std::vector<double> efields = m_eFields;
1855 for (auto& field : efields) field /= m_pressureTable;
1856 int cnt = 0;
1857 // List 5 values, then new line.
1858 PrintArray(efields, outfile, cnt, 5);
1859 if (nE % 5 != 0) outfile << "\n";
1860 // Store angles.
1861 outfile << " E-B angles \n";
1862 cnt = 0;
1863 PrintArray(m_bAngles, outfile, cnt, 5);
1864 if (nA % 5 != 0) outfile << "\n";
1865 // Store B fields (convert to hGauss).
1866 outfile << " B fields \n";
1867 std::vector<double> bfields = m_bFields;
1868 for (auto& field : bfields) field *= 100.;
1869 cnt = 0;
1870 PrintArray(bfields, outfile, cnt, 5);
1871 if (nB % 5 != 0) outfile << "\n";
1872
1873 // Store the gas composition.
1874 outfile << " Mixture: \n";
1875 cnt = 0;
1876 PrintArray(mixture, outfile, cnt, 5);
1877 if (nMagboltzGases % 5 != 0) outfile << "\n";
1878
1879 cnt = 0;
1880 for (const auto& exc : m_excLevels) {
1881 ++cnt;
1882 outfile << " Excitation " << FmtInt(cnt, 5) << ": " << std::setw(45);
1883 // If the label contains white space, enclose it in quotes.
1884 if (exc.label.find(" ") != std::string::npos) {
1885 outfile << "\"" + exc.label + "\"";
1886 } else {
1887 outfile << exc.label;
1888 }
1889 outfile << " " << FmtFloat(exc.energy) << FmtFloat(exc.prob)
1890 << FmtFloat(exc.rms) << FmtFloat(exc.dt) << "\n";
1891 }
1892 cnt = 0;
1893 for (const auto& ion : m_ionLevels) {
1894 ++cnt;
1895 outfile << " Ionisation " << FmtInt(cnt, 5) << ": " << std::setw(45);
1896 // If the label contains white space, enclose it in quotes.
1897 if (ion.label.find(" ") != std::string::npos) {
1898 outfile << "\"" + ion.label << "\"";
1899 } else {
1900 outfile << ion.label;
1901 }
1902 outfile << " " << FmtFloat(ion.energy) << "\n";
1903 }
1904
1905 const double sqrp = sqrt(m_pressureTable);
1906 const double logp = log(m_pressureTable);
1907 outfile << " The gas tables follow:\n";
1908 cnt = 0;
1909 for (unsigned int i = 0; i < nE; ++i) {
1910 for (unsigned int j = 0; j < nA; ++j) {
1911 for (unsigned int k = 0; k < nB; ++k) {
1912 // Get the velocities.
1913 double ve = m_eVelE.empty() ? 0. : m_eVelE[j][k][i];
1914 double vb = m_eVelB.empty() ? 0. : m_eVelB[j][k][i];
1915 double vx = m_eVelX.empty() ? 0. : m_eVelX[j][k][i];
1916 // Convert from cm / ns to cm / us.
1917 ve *= 1.e3;
1918 vb *= 1.e3;
1919 vx *= 1.e3;
1920 // Make a list of the values to be written, start with the velocities.
1921 std::vector<double> val;
1922 if (m_tab2d) {
1923 val = {ve, vb, vx};
1924 } else {
1925 // Add dummy spline values in case of a 1D table.
1926 val = {ve, 0., vb, 0., vx, 0.};
1927 }
1928 // Get the diffusion coefficients.
1929 double dl = m_eDifL.empty() ? 0. : m_eDifL[j][k][i] * sqrp;
1930 double dt = m_eDifT.empty() ? 0. : m_eDifT[j][k][i] * sqrp;
1931 // Get the Townsend and attachment coefficients.
1932 double alpha = m_eAlp.empty() ? -30. : m_eAlp[j][k][i] - logp;
1933 double alpha0 = m_eAlp0.empty() ? -30. : m_eAlp0[j][k][i] - logp;
1934 double eta = m_eAtt.empty() ? -30. : m_eAtt[j][k][i] - logp;
1935 // Add them to the list.
1936 if (m_tab2d) {
1937 val.insert(val.end(), {dl, dt, alpha, alpha0, eta});
1938 } else {
1939 val.insert(val.end(), {dl, 0., dt, 0., alpha, 0., alpha0, eta, 0.});
1940 }
1941 // Get the ion mobility and convert from cm2 / (V ns) to cm2 / (V us).
1942 double mu = m_iMob.empty() ? 0. : 1.e3 * m_iMob[j][k][i];
1943 // Get the Lorentz angle.
1944 double lor = m_eLor.empty() ? 0 : m_eLor[j][k][i];
1945 // Get the dissociation coefficient.
1946 double diss = m_iDis.empty() ? -30. : m_iDis[j][k][i] - logp;
1947 // Add them to the list.
1948 if (m_tab2d) {
1949 val.insert(val.end(), {mu, lor, diss});
1950 } else {
1951 val.insert(val.end(), {mu, 0., lor, 0., diss, 0.});
1952 }
1953 // Get the components of the diffusion tensor.
1954 for (int l = 0; l < 6; ++l) {
1955 if (!m_eDifM.empty()) {
1956 const double cov = m_eDifM[l][j][k][i] * m_pressureTable;
1957 val.push_back(cov);
1958 } else {
1959 val.push_back(0.);
1960 }
1961 if (!m_tab2d) val.push_back(0.);
1962 }
1963 // Get the excitation and ionisation rates.
1964 for (const auto& rexc : m_excRates) {
1965 if (rexc[j][k][i] > Small) {
1966 val.push_back(rexc[j][k][i]);
1967 } else {
1968 val.push_back(0.);
1969 }
1970 if (!m_tab2d) val.push_back(0.);
1971 }
1972 for (const auto& rion : m_ionRates) {
1973 if (rion[j][k][i] > Small) {
1974 val.push_back(rion[j][k][i]);
1975 } else {
1976 val.push_back(0.);
1977 }
1978 if (!m_tab2d) val.push_back(0.);
1979 }
1980 PrintArray(val, outfile, cnt, 8);
1981 }
1982 if (cnt % 8 != 0) outfile << "\n";
1983 cnt = 0;
1984 }
1985 }
1986
1987 if (!m_tab2d) {
1988 // Extrapolation methods
1989 int extrapH[13], extrapL[13];
1990 extrapL[0] = extrapL[1] = extrapL[2] = m_extrVel.first;
1991 extrapH[0] = extrapH[1] = extrapH[2] = m_extrVel.second;
1992 extrapL[3] = extrapL[8] = extrapL[10] = m_extrDif.first;
1993 extrapH[3] = extrapH[8] = extrapH[10] = m_extrDif.second;
1994 extrapL[4] = m_extrAlp.first;
1995 extrapH[4] = m_extrAlp.second;
1996 extrapL[5] = m_extrAtt.first;
1997 extrapH[5] = m_extrAtt.second;
1998 extrapL[6] = m_extrMob.first;
1999 extrapH[6] = m_extrMob.second;
2000 // Lorentz angle
2001 extrapL[7] = m_extrLor.first;
2002 extrapH[7] = m_extrLor.second;
2003 extrapL[9] = m_extrDis.first;
2004 extrapH[9] = m_extrDis.second;
2005 extrapL[11] = m_extrExc.first;
2006 extrapH[11] = m_extrExc.second;
2007 extrapL[12] = m_extrIon.first;
2008 extrapH[12] = m_extrIon.second;
2009 outfile << " H Extr: ";
2010 for (int i = 0; i < 13; i++) outfile << FmtInt(extrapH[i], 5);
2011 outfile << "\n";
2012 outfile << " L Extr: ";
2013 for (int i = 0; i < 13; i++) outfile << FmtInt(extrapL[i], 5);
2014 outfile << "\n";
2015 }
2016 // Increment the threshold indices for compatibility with Fortran.
2017 outfile << " Thresholds: " << FmtInt(m_eThrAlp + 1, 10)
2018 << FmtInt(m_eThrAtt + 1, 10) << FmtInt(m_iThrDis + 1, 10) << "\n";
2019 // Interpolation methods.
2020 int interp[13];
2021 interp[0] = interp[1] = interp[2] = m_intpVel;
2022 interp[3] = interp[8] = interp[10] = m_intpDif;
2023 interp[4] = m_intpAlp;
2024 interp[5] = m_intpAtt;
2025 interp[6] = m_intpMob;
2026 interp[7] = m_intpLor;
2027 interp[9] = m_intpDis;
2028 interp[11] = m_intpExc;
2029 interp[12] = m_intpIon;
2030 outfile << " Interp: ";
2031 for (int i = 0; i < 13; i++) outfile << FmtInt(interp[i], 5);
2032 outfile << "\n";
2033 outfile << " A =" << FmtFloat(0.) << ", Z =" << FmtFloat(0.) << ","
2034 << " EMPROB=" << FmtFloat(0.) << ", EPAIR =" << FmtFloat(0.) << "\n";
2035 const double dli = m_iDifL.empty() ? 0. : m_iDifL[0][0][0];
2036 const double dti = m_iDifT.empty() ? 0. : m_iDifT[0][0][0];
2037 outfile << " Ion diffusion: " << FmtFloat(dli) << FmtFloat(dti) << "\n";
2038 outfile << " CMEAN =" << FmtFloat(0.) << ", RHO =" << FmtFloat(0.) << ","
2039 << " PGAS =" << FmtFloat(m_pressureTable) << ","
2040 << " TGAS =" << FmtFloat(m_temperatureTable) << "\n";
2041 outfile << " CLSTYP : NOT SET \n"
2042 << " FCNCLS : " << std::string(80, ' ') << "\n"
2043 << " NCLS : " << FmtInt(0, 10) << "\n"
2044 << " Average : " << FmtFloat(0., 25, 18) << "\n"
2045 << " Heed initialisation done: F\n"
2046 << " SRIM initialisation done: F\n";
2047 outfile.close();
2048
2049 return true;
2050}
2051
2052void MediumGas::GetGasBits(std::bitset<20>& gasok) const {
2053
2054 gasok.reset();
2055 if (!m_eVelE.empty()) gasok.set(0);
2056 if (!m_iMob.empty()) gasok.set(1);
2057 if (!m_eDifL.empty()) gasok.set(2);
2058 if (!m_eAlp.empty()) gasok.set(3);
2059 // Cluster size distribution; skipped
2060 if (!m_eAtt.empty()) gasok.set(5);
2061 if (!m_eLor.empty()) gasok.set(6);
2062 if (!m_eDifT.empty()) gasok.set(7);
2063 if (!m_eVelB.empty()) gasok.set(8);
2064 if (!m_eVelX.empty()) gasok.set(9);
2065 if (!m_eDifM.empty()) gasok.set(10);
2066 if (!m_iDis.empty()) gasok.set(11);
2067 // SRIM, HEED; skipped
2068 if (!m_excRates.empty()) gasok.set(14);
2069 if (!m_ionRates.empty()) gasok.set(15);
2070}
2071
2073 // Print a summary.
2074 std::cout << m_className << "::PrintGas:\n"
2075 << " Gas composition: " << m_name;
2076 if (m_nComponents > 1) {
2077 std::cout << " (" << m_fraction[0] * 100;
2078 for (unsigned int i = 1; i < m_nComponents; ++i) {
2079 std::cout << "/" << m_fraction[i] * 100;
2080 }
2081 std::cout << ")";
2082 }
2083 std::cout << "\n";
2084 std::cout << " Pressure: " << m_pressure << " Torr\n"
2085 << " Temperature: " << m_temperature << " K\n"
2086 << " Gas file:\n"
2087 << " Pressure: " << m_pressureTable << " Torr\n"
2088 << " Temperature: " << m_temperatureTable << " K\n";
2089 if (m_eFields.size() > 1) {
2090 std::cout << " Electric field range: " << m_eFields[0] << " - "
2091 << m_eFields.back() << " V/cm in " << m_eFields.size() - 1
2092 << " steps.\n";
2093 } else if (m_eFields.size() == 1) {
2094 std::cout << " Electric field: " << m_eFields[0] << " V/cm\n";
2095 } else {
2096 std::cout << " Electric field range: not set\n";
2097 }
2098 if (m_bFields.size() > 1) {
2099 std::cout << " Magnetic field range: " << m_bFields[0] << " - "
2100 << m_bFields.back() << " T in " << m_bFields.size() - 1
2101 << " steps.\n";
2102 } else if (m_bFields.size() == 1) {
2103 std::cout << " Magnetic field: " << m_bFields[0] << " T\n";
2104 } else {
2105 std::cout << " Magnetic field range: not set\n";
2106 }
2107 if (m_bAngles.size() > 1) {
2108 std::cout << " Angular range: " << m_bAngles[0] << " - "
2109 << m_bAngles.back() << " rad in " << m_bAngles.size() - 1
2110 << " steps.\n";
2111 } else if (m_bAngles.size() == 1) {
2112 std::cout << " Angle between E and B: " << m_bAngles[0] << " rad\n";
2113 } else {
2114 std::cout << " Angular range: not set\n";
2115 }
2116
2117 std::cout << " Available electron transport data:\n";
2118 if (!m_eVelE.empty()) {
2119 std::cout << " Velocity along E\n";
2120 }
2121 if (!m_eVelB.empty()) {
2122 std::cout << " Velocity along Bt\n";
2123 }
2124 if (!m_eVelX.empty()) {
2125 std::cout << " Velocity along ExB\n";
2126 }
2127 if (!(m_eVelE.empty() && m_eVelB.empty() &&
2128 m_eVelX.empty())) {
2129 PrintExtrapolation(m_extrVel);
2130 std::cout << " Interpolation order: " << m_intpVel << "\n";
2131 }
2132 if (!m_eDifL.empty()) {
2133 std::cout << " Longitudinal diffusion coefficient\n";
2134 }
2135 if (!m_eDifT.empty()) {
2136 std::cout << " Transverse diffusion coefficient\n";
2137 }
2138 if (!m_eDifM.empty()) {
2139 std::cout << " Diffusion tensor\n";
2140 }
2141 if (!m_eDifL.empty() || !m_eDifT.empty() || !m_eDifM.empty()) {
2142 PrintExtrapolation(m_extrDif);
2143 std::cout << " Interpolation order: " << m_intpDif << "\n";
2144 }
2145 if (!m_eAlp.empty()) {
2146 std::cout << " Townsend coefficient\n";
2147 PrintExtrapolation(m_extrAlp);
2148 std::cout << " Interpolation order: " << m_intpAlp << "\n";
2149 }
2150 if (!m_eAtt.empty()) {
2151 std::cout << " Attachment coefficient\n";
2152 PrintExtrapolation(m_extrAtt);
2153 std::cout << " Interpolation order: " << m_intpAtt << "\n";
2154 }
2155 if (!m_eLor.empty()) {
2156 std::cout << " Lorentz Angle\n";
2157 PrintExtrapolation(m_extrLor);
2158 std::cout << " Interpolation order: " << m_intpLor << "\n";
2159 }
2160 if (!m_excRates.empty()) {
2161 std::cout << " Excitation rates\n";
2162 for (const auto& exc : m_excLevels) {
2163 std::cout << " " << exc.label << "\n";
2164 std::cout << " Energy = " << exc.energy << " eV";
2165 if (exc.prob > 0.) {
2166 std::cout << ", Penning transfer probability = " << exc.prob;
2167 }
2168 std::cout << "\n";
2169 }
2170 PrintExtrapolation(m_extrExc);
2171 std::cout << " Interpolation order: " << m_intpExc << "\n";
2172 }
2173 if (!m_ionRates.empty()) {
2174 std::cout << " Ionisation rates\n";
2175 for (const auto& ion : m_ionLevels) {
2176 std::cout << " " << ion.label << "\n";
2177 std::cout << " Threshold = " << ion.energy << " eV\n";
2178 }
2179 PrintExtrapolation(m_extrIon);
2180 std::cout << " Interpolation order: " << m_intpIon << "\n";
2181 }
2182 if (m_eVelE.empty() && m_eVelB.empty() && m_eVelX.empty() &&
2183 m_eDifL.empty() && m_eDifT.empty() && m_eDifM.empty() &&
2184 m_eAlp.empty() && m_eAtt.empty() && m_excRates.empty() &&
2185 m_ionRates.empty() && m_eLor.empty()) {
2186 std::cout << " none\n";
2187 }
2188
2189 std::cout << " Available ion transport data:\n";
2190 if (!m_iMob.empty()) {
2191 std::cout << " Mobility\n";
2192 PrintExtrapolation(m_extrMob);
2193 std::cout << " Interpolation order: " << m_intpMob << "\n";
2194 }
2195 if (!m_iDifL.empty()) {
2196 std::cout << " Longitudinal diffusion coefficient\n";
2197 }
2198 if (!m_iDifT.empty()) {
2199 std::cout << " Transverse diffusion coefficient\n";
2200 }
2201 if (!m_iDifL.empty() || !m_iDifT.empty()) {
2202 PrintExtrapolation(m_extrDif);
2203 std::cout << " Interpolation order: " << m_intpDif << "\n";
2204 }
2205 if (!m_iDis.empty()) {
2206 std::cout << " Dissociation coefficient\n";
2207 PrintExtrapolation(m_extrDis);
2208 std::cout << " Interpolation order: " << m_intpDis << "\n";
2209 }
2210 if (m_iMob.empty() && m_iDifL.empty() && m_iDifT.empty() && m_iDis.empty()) {
2211 std::cout << " none\n";
2212 }
2213}
2214
2215bool MediumGas::LoadIonMobility(const std::string& filename) {
2216 // Open the file.
2217 std::ifstream infile;
2218 infile.open(filename.c_str(), std::ios::in);
2219 // Make sure the file could actually be opened.
2220 if (!infile) {
2221 std::cerr << m_className << "::LoadIonMobility:\n"
2222 << " Error opening file " << filename << ".\n";
2223 return false;
2224 }
2225
2226 std::vector<std::pair<double, double> > data;
2227 // Read the file line by line.
2228 int i = 0;
2229 constexpr size_t size = 100;
2230 char line[size];
2231 while (infile.getline(line, size)) {
2232 ++i;
2233 char* token = strtok(line, " ,\t");
2234 if (!token) break;
2235 if (strcmp(token, "#") == 0 || strcmp(token, "*") == 0 ||
2236 strcmp(token, "//") == 0) {
2237 continue;
2238 }
2239 double field = atof(token);
2240 token = strtok(NULL, " ,\t");
2241 if (!token) {
2242 std::cerr << m_className << "::LoadIonMobility:\n"
2243 << " Found E/N but no mobility before the end-of-line.\n"
2244 << " Skipping line " << i << ".\n";
2245 continue;
2246 }
2247 double mu = atof(token);
2248 if (m_debug) {
2249 std::cout << " E/N = " << field << " Td: mu = " << mu << " cm2/(Vs)\n";
2250 }
2251 // Make sure the values make sense.
2252 // Negative field values are not allowed.
2253 if (field < 0.) {
2254 std::cerr << m_className << "::LoadIonMobility:\n"
2255 << " Negative electric field (line " << i << ").\n";
2256 return false;
2257 }
2258 // Add the values to the list.
2259 data.push_back(std::make_pair(field, mu));
2260 }
2261 infile.close();
2262
2263 if (data.empty()) {
2264 std::cerr << m_className << "::LoadIonMobilities:\n"
2265 << " No valid data found.\n";
2266 return false;
2267 }
2268 // Sort by electric field.
2269 std::sort(data.begin(), data.end());
2270
2271 // The E/N values in the file are supposed to be in Td (10^-17 V cm2).
2272 const double scaleField = 1.e-17 * GetNumberDensity();
2273 // The reduced mobilities in the file are supposed to be in cm2/(V s).
2274 const double scaleMobility = 1.e-9 * (AtmosphericPressure / m_pressure) *
2275 (m_temperature / ZeroCelsius);
2276
2277 const size_t ne = data.size();
2278 std::vector<double> efields(ne, 0.);
2279 std::vector<double> mobilities(ne, 0.);
2280 for (size_t j = 0; j < ne; ++j) {
2281 // Scale the fields and mobilities.
2282 efields[j] = data[j].first * scaleField;
2283 mobilities[j] = data[j].second * scaleMobility;
2284 }
2285
2286 std::cout << m_className << "::LoadIonMobility:\n"
2287 << " Read " << ne << " values from file " << filename << "\n";
2288
2289 return SetIonMobility(efields, mobilities);
2290}
2291
2293
2295 m_eAlp0.clear();
2296 m_excLevels.clear();
2297 m_ionLevels.clear();
2298 m_excRates.clear();
2299 m_ionRates.clear();
2300}
2301
2303 const double lambda) {
2304
2305 if (r < 0. || r > 1.) {
2306 std::cerr << m_className << "::EnablePenningTransfer:\n"
2307 << " Transfer probability must be in the range [0, 1].\n";
2308 return false;
2309 }
2310
2311 m_rPenningGlobal = r;
2312 m_lambdaPenningGlobal = lambda > Small ? lambda : 0.;
2313
2314 std::cout << m_className << "::EnablePenningTransfer:\n"
2315 << " Global Penning transfer parameters set to:\n"
2316 << " r = " << m_rPenningGlobal << "\n"
2317 << " lambda = " << m_lambdaPenningGlobal << " cm\n";
2318
2319 // Find the min. ionisation energy.
2320 if (m_ionLevels.empty()) {
2321 std::cerr << m_className << "::EnablePenningTransfer:\n Warning: present"
2322 << " gas table has no ionisation rates.\n Ignore this message "
2323 << "if you are using microscopic tracking only.\n";
2324 return true;
2325 }
2326 double minIonPot = -1.;
2327 for (const auto& ion : m_ionLevels) {
2328 if (minIonPot < 0.) {
2329 minIonPot = ion.energy;
2330 } else {
2331 minIonPot = std::min(minIonPot, ion.energy);
2332 }
2333 }
2334
2335 // Update the transfer probabilities of the excitation levels in the table.
2336 unsigned int nLevelsFound = 0;
2337 for (auto& exc : m_excLevels) {
2338 if (exc.energy < minIonPot) continue;
2339 exc.prob = m_rPenningGlobal;
2340 exc.rms = m_lambdaPenningGlobal;
2341 ++nLevelsFound;
2342 }
2343 if (nLevelsFound > 0) {
2344 std::cout << m_className << "::EnablePenningTransfer:\n"
2345 << " Updated transfer probabilities for " << nLevelsFound
2346 << " excitation rates.\n";
2348 } else {
2349 std::cerr << m_className << "::EnablePenningTransfer:\n Warning: present"
2350 << " gas table has no eligible excitation rates.\n Ignore this"
2351 << " message if you are using microscopic tracking only.\n";
2352 }
2353 return true;
2354}
2355
2356bool MediumGas::EnablePenningTransfer(const double r, const double lambda,
2357 std::string gasname) {
2358
2359 if (r < 0. || r > 1.) {
2360 std::cerr << m_className << "::EnablePenningTransfer:\n"
2361 << " Transfer probability must be in the range [0, 1].\n";
2362 return false;
2363 }
2364
2365 // Get the "standard" name of this gas.
2366 gasname = GetGasName(gasname);
2367 if (gasname.empty()) {
2368 std::cerr << m_className << "::EnablePenningTransfer: Unknown gas name.\n";
2369 return false;
2370 }
2371
2372 // Look for this gas in the present gas mixture.
2373 int iGas = -1;
2374 for (unsigned int i = 0; i < m_nComponents; ++i) {
2375 if (m_gas[i] == gasname) {
2376 m_rPenningGas[i] = r;
2377 m_lambdaPenningGas[i] = lambda > Small ? lambda : 0.;
2378 iGas = i;
2379 break;
2380 }
2381 }
2382
2383 if (iGas < 0) {
2384 std::cerr << m_className << "::EnablePenningTransfer:\n"
2385 << " Requested gas (" << gasname
2386 << ") is not part of the present gas mixture.\n";
2387 return false;
2388 }
2389
2390 // Find the min. ionisation energy.
2391 if (m_ionLevels.empty()) {
2392 std::cerr << m_className << "::EnablePenningTransfer:\n Warning: present"
2393 << " gas table has no ionisation rates.\n Ignore this message"
2394 << " if you are using microscopic tracking only.\n";
2395 return true;
2396 }
2397 double minIonPot = -1.;
2398 for (const auto& ion : m_ionLevels) {
2399 if (minIonPot < 0.) {
2400 minIonPot = ion.energy;
2401 } else {
2402 minIonPot = std::min(minIonPot, ion.energy);
2403 }
2404 }
2405 // Update the transfer probabilities of the excitation levels in the table.
2406 unsigned int nLevelsFound = 0;
2407 for (auto& exc : m_excLevels) {
2408 if (exc.energy < minIonPot) continue;
2409 // Try to extract the gas name from the label.
2410 // TODO: test if this works for all gases and excitation levels.
2411 const auto pos = exc.label.find('-');
2412 if (pos == std::string::npos) continue;
2413 if (GetGasName(exc.label.substr(0, pos)) != gasname) continue;
2414 exc.prob = r;
2415 exc.rms = lambda;
2416 ++nLevelsFound;
2417 }
2418 if (nLevelsFound > 0) {
2419 std::cout << m_className << "::EnablePenningTransfer:\n"
2420 << " Updated transfer probabilities for " << nLevelsFound
2421 << " excitation rates.\n";
2423 } else {
2424 std::cerr << m_className << "::EnablePenningTransfer:\n Warning: present"
2425 << " gas table has no eligible excitation rates.\n Ignore this"
2426 << " message if you are using microscopic tracking only.\n";
2427 }
2428 return true;
2429}
2430
2432
2433 m_rPenningGlobal = 0.;
2435
2436 m_rPenningGas.fill(0.);
2437 m_lambdaPenningGas.fill(0.);
2438
2439 if (m_excLevels.empty()) return;
2440 for (auto& exc : m_excLevels) {
2441 exc.prob = 0.;
2442 }
2444}
2445
2446bool MediumGas::DisablePenningTransfer(std::string gasname) {
2447
2448 // Get the "standard" name of this gas.
2449 gasname = GetGasName(gasname);
2450 if (gasname.empty()) {
2451 std::cerr << m_className << "::DisablePenningTransfer: Unknown gas name.\n";
2452 return false;
2453 }
2454
2455 // Look for this gas in the present gas mixture.
2456 int iGas = -1;
2457 for (unsigned int i = 0; i < m_nComponents; ++i) {
2458 if (m_gas[i] == gasname) {
2459 m_rPenningGas[i] = 0.;
2460 m_lambdaPenningGas[i] = 0.;
2461 iGas = i;
2462 break;
2463 }
2464 }
2465
2466 if (iGas < 0) {
2467 std::cerr << m_className << "::DisablePenningTransfer:\n"
2468 << " Requested gas (" << gasname
2469 << ") is not part of the present gas mixture.\n";
2470 return false;
2471 }
2472
2473 if (m_excLevels.empty()) return true;
2474 for (auto& exc : m_excLevels) {
2475 // Try to extract the gas name from the label.
2476 const auto pos = exc.label.find('-');
2477 if (pos == std::string::npos) continue;
2478 if (GetGasName(exc.label.substr(0, pos)) != gasname) continue;
2479 exc.prob = 0.;
2480 }
2482 return true;
2483}
2484
2485
2487
2488 // -----------------------------------------------------------------------
2489 // GASSPT
2490 // -----------------------------------------------------------------------
2491
2492 // Make sure there are Townsend coefficients.
2493 if (m_eAlp.empty() || m_eAlp0.empty()) {
2494 std::cerr << m_className << "::AdjustTownsendCoefficient:\n "
2495 << "Present gas table does not include Townsend coefficients.\n";
2496 return false;
2497 }
2498 // Make sure there are excitation and ionisation rates.
2499 if (m_excLevels.empty() || m_excRates.empty()) {
2500 std::cerr << m_className << "::AdjustTownsendCoefficient:\n "
2501 << "Present gas table does not include excitation rates.\n";
2502 return false;
2503 }
2504 if (m_ionLevels.empty() || m_ionRates.empty()) {
2505 std::cerr << m_className << "::AdjustTownsendCoefficient:\n "
2506 << "Present gas table does not include ionisation rates.\n";
2507 return false;
2508 }
2509 const unsigned int nE = m_eFields.size();
2510 const unsigned int nB = m_bFields.size();
2511 const unsigned int nA = m_bAngles.size();
2512 if (m_debug) {
2513 std::cout << m_className << "::AdjustTownsendCoefficient:\n"
2514 << " Entry Exc. Ion.\n";
2515 }
2516 for (unsigned int i = 0; i < nE; ++i) {
2517 for (unsigned int j = 0; j < nA; ++j) {
2518 for (unsigned int k = 0; k < nB; ++k) {
2519 // Compute total ionisation rate.
2520 double rion = 0.;
2521 for (const auto& ion : m_ionRates) {
2522 rion += ion[j][k][i];
2523 }
2524 // Compute rate of Penning ionisations.
2525 double rexc = 0.;
2526 const unsigned int nexc = m_excLevels.size();
2527 for (unsigned int ie = 0; ie < nexc; ++ie) {
2528 rexc += m_excLevels[ie].prob * m_excRates[ie][j][k][i];
2529 }
2530 if (m_debug) {
2531 std::cout << FmtInt(i, 4) << FmtInt(j, 4) << FmtInt(k, 4)
2532 << FmtFloat(rexc, 12, 5) << FmtFloat(rion, 12, 5) << "\n";
2533 }
2534 // Adjust the Townsend coefficient.
2535 double alpha0 = m_eAlp0[j][k][i];
2536 if (alpha0 < -20.) {
2537 alpha0 = 0.;
2538 } else {
2539 alpha0 = m_pressure * exp(alpha0);
2540 }
2541 double alpha1 = alpha0;
2542 if (rion > 0.) alpha1 *= (rexc + rion) / rion;
2543 m_eAlp[j][k][i] = alpha1 > 0. ? log(alpha1 / m_pressure) : -30.;
2544 }
2545 }
2546 }
2547 // Update the threshold index.
2549 return true;
2550}
2551
2552bool MediumGas::GetGasInfo(const std::string& gasname, double& a,
2553 double& z) const {
2554 if (gasname == "CF4") {
2555 a = 12.0107 + 4 * 18.9984032;
2556 z = 6 + 4 * 9;
2557 return true;
2558 } else if (gasname == "Ar") {
2559 a = 39.948;
2560 z = 18;
2561 } else if (gasname == "He") {
2562 a = 4.002602;
2563 z = 2;
2564 } else if (gasname == "He-3") {
2565 a = 3.01602931914;
2566 z = 2;
2567 } else if (gasname == "Ne") {
2568 a = 20.1797;
2569 z = 10;
2570 } else if (gasname == "Kr") {
2571 a = 37.798;
2572 z = 36;
2573 } else if (gasname == "Xe") {
2574 a = 131.293;
2575 z = 54;
2576 } else if (gasname == "CH4") {
2577 a = 12.0107 + 4 * 1.00794;
2578 z = 6 + 4;
2579 } else if (gasname == "C2H6") {
2580 a = 2 * 12.0107 + 6 * 1.00794;
2581 z = 2 * 6 + 6;
2582 } else if (gasname == "C3H8") {
2583 a = 3 * 12.0107 + 8 * 1.00794;
2584 z = 3 * 6 + 8;
2585 } else if (gasname == "iC4H10") {
2586 a = 4 * 12.0107 + 10 * 1.00794;
2587 z = 4 * 6 + 10;
2588 } else if (gasname == "CO2") {
2589 a = 12.0107 + 2 * 15.9994;
2590 z = 6 + 2 * 8;
2591 } else if (gasname == "neoC5H12") {
2592 a = 5 * 12.0107 + 12 * 1.00794;
2593 z = 5 * 6 + 12;
2594 } else if (gasname == "H2O") {
2595 a = 2 * 1.00794 + 15.9994;
2596 z = 2 + 8;
2597 } else if (gasname == "O2") {
2598 a = 2 * 15.9994;
2599 z = 2 * 8;
2600 } else if (gasname == "N2") {
2601 a = 2 * 14.0067;
2602 z = 2 * 7;
2603 } else if (gasname == "NO") {
2604 a = 14.0067 + 15.9994;
2605 z = 7 + 8;
2606 } else if (gasname == "N2O") {
2607 a = 2 * 14.0067 + 15.9994;
2608 z = 2 * 7 + 8;
2609 } else if (gasname == "C2H4") {
2610 a = 2 * 12.0107 + 4 * 1.00794;
2611 z = 2 * 6 + 4;
2612 } else if (gasname == "C2H2") {
2613 a = 2 * 12.0107 + 2 * 1.00794;
2614 z = 2 * 6 + 2;
2615 } else if (gasname == "H2" || gasname == "paraH2") {
2616 a = 2 * 1.00794;
2617 z = 2;
2618 } else if (gasname == "D2" || gasname == "orthoD2") {
2619 a = 2 * 2.01410177785;
2620 z = 2;
2621 } else if (gasname == "CO") {
2622 a = 12.0107 + 15.9994;
2623 z = 6 + 8;
2624 } else if (gasname == "Methylal") {
2625 a = 3 * 12.0107 + 8 * 1.00794 + 2 * 15.9994;
2626 z = 3 * 6 + 8 + 2 * 8;
2627 } else if (gasname == "DME") {
2628 a = 4 * 12.0107 + 10 * 1.00794 + 2 * 15.9994;
2629 z = 4 * 6 + 10 + 2 * 8;
2630 } else if (gasname == "Reid-Step" || gasname == "Maxwell-Model" ||
2631 gasname == "Reid-Ramp") {
2632 a = 1.;
2633 z = 1.;
2634 } else if (gasname == "C2F6") {
2635 a = 2 * 12.0107 + 6 * 18.9984032;
2636 z = 2 * 6 + 6 * 9;
2637 } else if (gasname == "SF6") {
2638 a = 32.065 + 6 * 18.9984032;
2639 z = 16 + 6 * 9;
2640 } else if (gasname == "NH3") {
2641 a = 14.0067 + 3 * 1.00794;
2642 z = 7 + 3;
2643 } else if (gasname == "C3H6") {
2644 a = 3 * 12.0107 + 6 * 1.00794;
2645 z = 3 * 6 + 6;
2646 } else if (gasname == "cC3H6") {
2647 a = 3 * 12.0107 + 6 * 1.00794;
2648 z = 3 * 6 + 6;
2649 } else if (gasname == "CH3OH") {
2650 a = 12.0107 + 4 * 1.00794 + 15.9994;
2651 z = 6 + 4 + 8;
2652 } else if (gasname == "C2H5OH") {
2653 a = 2 * 12.0107 + 6 * 1.00794 + 15.9994;
2654 z = 2 * 6 + 6 + 8;
2655 } else if (gasname == "C3H7OH" || gasname == "nC3H7OH") {
2656 a = 3 * 12.0107 + 8 * 1.00794 + 15.9994;
2657 z = 3 * 6 + 8 * 8;
2658 } else if (gasname == "Cs") {
2659 a = 132.9054519;
2660 z = 55;
2661 } else if (gasname == "F2") {
2662 a = 2 * 18.9984032;
2663 z = 2 * 9;
2664 } else if (gasname == "CS2") {
2665 a = 12.0107 + 2 * 32.065;
2666 z = 6 + 2 * 16;
2667 } else if (gasname == "COS") {
2668 a = 12.0107 + 15.9994 + 32.065;
2669 z = 6 + 8 + 16;
2670 } else if (gasname == "CD4") {
2671 a = 12.0107 + 4 * 2.01410177785;
2672 z = 6 + 4;
2673 } else if (gasname == "BF3") {
2674 a = 10.811 + 3 * 18.9984032;
2675 z = 5 + 3 * 9;
2676 } else if (gasname == "C2H2F4") {
2677 a = 2 * 12.0107 + 2 * 1.00794 + 4 * 18.9984032;
2678 z = 2 * 6 + 2 + 4 * 9;
2679 } else if (gasname == "CHF3") {
2680 a = 12.0107 + 1.00794 + 3 * 18.9984032;
2681 z = 6 + 1 + 3 * 9;
2682 } else if (gasname == "CF3Br") {
2683 a = 12.0107 + 3 * 18.9984032 + 79.904;
2684 z = 6 + 3 * 9 + 35;
2685 } else if (gasname == "C3F8") {
2686 a = 3 * 12.0107 + 8 * 18.9984032;
2687 z = 3 * 6 + 8 * 9;
2688 } else if (gasname == "O3") {
2689 a = 3 * 15.9994;
2690 z = 3 * 8;
2691 } else if (gasname == "Hg") {
2692 a = 2 * 200.59;
2693 z = 80;
2694 } else if (gasname == "H2S") {
2695 a = 2 * 1.00794 + 32.065;
2696 z = 2 + 16;
2697 } else if (gasname == "nC4H10") {
2698 a = 4 * 12.0107 + 10 * 1.00794;
2699 z = 4 * 6 + 10;
2700 } else if (gasname == "nC5H12") {
2701 a = 5 * 12.0107 + 12 * 1.00794;
2702 z = 5 * 6 + 12;
2703 } else if (gasname == "GeH4") {
2704 a = 72.64 + 4 * 1.00794;
2705 z = 32 + 4;
2706 } else if (gasname == "SiH4") {
2707 a = 28.0855 + 4 * 1.00794;
2708 z = 14 + 4;
2709 } else {
2710 a = 0.;
2711 z = 0.;
2712 return false;
2713 }
2714
2715 return true;
2716}
2717
2718std::string MediumGas::GetGasName(const int gasnumber, const int version) const {
2719
2720 switch (gasnumber) {
2721 case 1:
2722 return "CF4";
2723 case 2:
2724 return "Ar";
2725 case 3:
2726 return "He";
2727 case 4:
2728 return "He-3";
2729 case 5:
2730 return "Ne";
2731 case 6:
2732 return "Kr";
2733 case 7:
2734 return "Xe";
2735 case 8:
2736 return "CH4";
2737 case 9:
2738 return "C2H6";
2739 case 10:
2740 return "C3H8";
2741 case 11:
2742 return "iC4H10";
2743 case 12:
2744 return "CO2";
2745 case 13:
2746 return "neoC5H12";
2747 case 14:
2748 return "H2O";
2749 case 15:
2750 return "O2";
2751 case 16:
2752 return "N2";
2753 case 17:
2754 return "NO";
2755 case 18:
2756 return "N2O";
2757 case 19:
2758 return "C2H4";
2759 case 20:
2760 return "C2H2";
2761 case 21:
2762 return "H2";
2763 case 22:
2764 return "D2";
2765 case 23:
2766 return "CO";
2767 case 24:
2768 return "Methylal";
2769 case 25:
2770 return "DME";
2771 case 26:
2772 return "Reid-Step";
2773 case 27:
2774 return "Maxwell-Model";
2775 case 28:
2776 return "Reid-Ramp";
2777 case 29:
2778 return "C2F6";
2779 case 30:
2780 return "SF6";
2781 case 31:
2782 return "NH3";
2783 case 32:
2784 return "C3H6";
2785 case 33:
2786 return "cC3H6";
2787 case 34:
2788 return "CH3OH";
2789 case 35:
2790 return "C2H5OH";
2791 case 36:
2792 return "C3H7OH";
2793 case 37:
2794 return "Cs";
2795 case 38:
2796 return "F2";
2797 case 39:
2798 return "CS2";
2799 case 40:
2800 return "COS";
2801 case 41:
2802 return "CD4";
2803 case 42:
2804 return "BF3";
2805 case 43:
2806 return "C2H2F4";
2807 case 44:
2808 return version <= 11 ? "He-3" : "TMA";
2809 case 45:
2810 return version <= 11 ? "He" : "paraH2";
2811 case 46:
2812 return version <= 11 ? "Ne" : "nC3H7OH";
2813 case 47:
2814 return "Ar";
2815 case 48:
2816 return version <= 11 ? "Kr" : "orthoD2";
2817 case 49:
2818 return "Xe";
2819 case 50:
2820 return "CHF3";
2821 case 51:
2822 return "CF3Br";
2823 case 52:
2824 return "C3F8";
2825 case 53:
2826 return "O3";
2827 case 54:
2828 return "Hg";
2829 case 55:
2830 return "H2S";
2831 case 56:
2832 return "nC4H10";
2833 case 57:
2834 return "nC5H12";
2835 case 58:
2836 return "N2";
2837 case 59:
2838 return "GeH4";
2839 case 60:
2840 return "SiH4";
2841 default:
2842 break;
2843 }
2844 return "";
2845}
2846
2847std::string MediumGas::GetGasName(std::string input) const {
2848 // Convert to upper-case.
2849 std::transform(input.begin(), input.end(), input.begin(), toupper);
2850
2851 if (input.empty()) return "";
2852
2853 if (input == "CF4" || input == "FREON" || input == "FREON-14" ||
2854 input == "TETRAFLUOROMETHANE") {
2855 return "CF4";
2856 } else if (input == "AR" || input == "ARGON") {
2857 return "Ar";
2858 } else if (input == "HE" || input == "HELIUM" || input == "HE-4" ||
2859 input == "HE 4" || input == "HE4" || input == "4-HE" ||
2860 input == "4 HE" || input == "4HE" || input == "HELIUM-4" ||
2861 input == "HELIUM 4" || input == "HELIUM4") {
2862 return "He";
2863 } else if (input == "HE-3" || input == "HE3" || input == "HELIUM-3" ||
2864 input == "HELIUM 3" || input == "HELIUM3") {
2865 return "He-3";
2866 } else if (input == "NE" || input == "NEON") {
2867 return "Ne";
2868 } else if (input == "KR" || input == "KRYPTON") {
2869 return "Kr";
2870 } else if (input == "XE" || input == "XENON") {
2871 return "Xe";
2872 } else if (input == "CH4" || input == "METHANE") {
2873 return "CH4";
2874 } else if (input == "C2H6" || input == "ETHANE") {
2875 return "C2H6";
2876 } else if (input == "C3H8" || input == "PROPANE") {
2877 return "C3H8";
2878 } else if (input == "C4H10" || input == "ISOBUTANE" || input == "ISO" ||
2879 input == "IC4H10" || input == "ISO-C4H10" || input == "ISOC4H10") {
2880 return "iC4H10";
2881 } else if (input == "CO2" || input == "CARBON-DIOXIDE" ||
2882 input == "CARBON DIOXIDE" || input == "CARBONDIOXIDE") {
2883 return "CO2";
2884 } else if (input == "NEOPENTANE" || input == "NEO-PENTANE" ||
2885 input == "NEO-C5H12" || input == "NEOC5H12" ||
2886 input == "DIMETHYLPROPANE" || input == "C5H12") {
2887 return "neoC5H12";
2888 } else if (input == "H2O" || input == "WATER" || input == "WATER-VAPOUR" ||
2889 input == "WATER VAPOUR") {
2890 return "H2O";
2891 } else if (input == "O2" || input == "OXYGEN") {
2892 return "O2";
2893 } else if (input == "NI" || input == "NITRO" || input == "N2" ||
2894 input == "NITROGEN") {
2895 return "N2";
2896 } else if (input == "NO" || input == "NITRIC-OXIDE" || input == "NITRIC OXIDE" ||
2897 input == "NITROGEN-MONOXIDE" || input == "NITROGEN MONOXIDE") {
2898 return "NO";
2899 } else if (input == "N2O" || input == "NITROUS-OXIDE" || input == "NITROUS OXIDE" ||
2900 input == "DINITROGEN-MONOXIDE" || input == "LAUGHING-GAS") {
2901 return "N2O";
2902 } else if (input == "C2H4" || input == "ETHENE" || input == "ETHYLENE") {
2903 return "C2H4";
2904 } else if (input == "C2H2" || input == "ACETYL" || input == "ACETYLENE" ||
2905 input == "ETHYNE") {
2906 return "C2H2";
2907 } else if (input == "H2" || input == "HYDROGEN") {
2908 return "H2";
2909 } else if (input == "PARA H2" || input == "PARA-H2" ||
2910 input == "PARAH2" || input == "PARA HYDROGEN" ||
2911 input == "PARA-HYDROGEN" || input == "PARAHYDROGEN") {
2912 return "paraH2";
2913 } else if (input == "D2" || input == "DEUTERIUM") {
2914 return "D2";
2915 } else if (input == "ORTHO D2" || input == "ORTHO-D2" ||
2916 input == "ORTHOD2" || input == "ORTHO DEUTERIUM" ||
2917 input == "ORTHO-DEUTERIUM" || input == "ORTHODEUTERIUM") {
2918 return "orthoD2";
2919 } else if (input == "CO" || input == "CARBON-MONOXIDE" ||
2920 input == "CARBON MONOXIDE") {
2921 return "CO";
2922 } else if (input == "METHYLAL" || input == "METHYLAL-HOT" || input == "DMM" ||
2923 input == "DIMETHOXYMETHANE" || input == "FORMAL" || input == "C3H8O2") {
2924 // Methylal (dimethoxymethane, CH3-O-CH2-O-CH3, "hot" version)
2925 return "Methylal";
2926 } else if (input == "DME" || input == "DIMETHYL-ETHER" || input == "DIMETHYLETHER" ||
2927 input == "DIMETHYL ETHER" || input == "METHYL ETHER" ||
2928 input == "METHYL-ETHER" || input == "METHYLETHER" ||
2929 input == "WOOD-ETHER" || input == "WOODETHER" || input == "WOOD ETHER" ||
2930 input == "DIMETHYL OXIDE" || input == "DIMETHYL-OXIDE" ||
2931 input == "DEMEON" || input == "METHOXYMETHANE" || input == "C4H10O2") {
2932 return "DME";
2933 } else if (input == "REID-STEP") {
2934 return "Reid-Step";
2935 } else if (input == "MAXWELL-MODEL") {
2936 return "Maxwell-Model";
2937 } else if (input == "REID-RAMP") {
2938 return "Reid-Ramp";
2939 } else if (input == "C2F6" || input == "FREON-116" || input == "ZYRON-116" ||
2940 input == "ZYRON-116-N5" || input == "HEXAFLUOROETHANE") {
2941 return "C2F6";
2942 } else if (input == "SF6" || input == "SULPHUR-HEXAFLUORIDE" ||
2943 input == "SULFUR-HEXAFLUORIDE" || input == "SULPHUR HEXAFLUORIDE" ||
2944 input == "SULFUR HEXAFLUORIDE") {
2945 return "SF6";
2946 } else if (input == "NH3" || input == "AMMONIA") {
2947 return "NH3";
2948 } else if (input == "C3H6" || input == "PROPENE" || input == "PROPYLENE") {
2949 return "C3H6";
2950 } else if (input == "C-PROPANE" || input == "CYCLO-PROPANE" ||
2951 input == "CYCLO PROPANE" || input == "CYCLOPROPANE" ||
2952 input == "C-C3H6" || input == "CC3H6" || input == "CYCLO-C3H6") {
2953 return "cC3H6";
2954 } else if (input == "METHANOL" || input == "METHYL-ALCOHOL" ||
2955 input == "METHYL ALCOHOL" || input == "WOOD ALCOHOL" ||
2956 input == "WOOD-ALCOHOL" || input == "CH3OH") {
2957 return "CH3OH";
2958 } else if (input == "ETHANOL" || input == "ETHYL-ALCOHOL" ||
2959 input == "ETHYL ALCOHOL" || input == "GRAIN ALCOHOL" ||
2960 input == "GRAIN-ALCOHOL" || input == "C2H5OH") {
2961 return "C2H5OH";
2962 } else if (input == "PROPANOL" || input == "2-PROPANOL" || input == "ISOPROPYL" ||
2963 input == "ISO-PROPANOL" || input == "ISOPROPANOL" ||
2964 input == "ISOPROPYL ALCOHOL" || input == "ISOPROPYL-ALCOHOL" ||
2965 input == "C3H7OH") {
2966 return "C3H7OH";
2967 } else if (input == "NPROPANOL" || input == "N-PROPANOL" ||
2968 input == "1-PROPANOL" || input == "PROPYL ALCOHOL" ||
2969 input == "PROPYL-ALCOHOL" || input == "N-PROPYL ALCOHOL" ||
2970 input == "NC3H7OH" || input == "N-C3H7OH") {
2971 return "nC3H7OH";
2972 } else if (input == "CS" || input == "CESIUM" || input == "CAESIUM") {
2973 return "Cs";
2974 } else if (input == "F2" || input == "FLUOR" || input == "FLUORINE") {
2975 return "F2";
2976 } else if (input == "CS2" || input == "CARBON-DISULPHIDE" ||
2977 input == "CARBON-DISULFIDE" || input == "CARBON DISULPHIDE" ||
2978 input == "CARBON DISULFIDE") {
2979 return "CS2";
2980 } else if (input == "COS" || input == "CARBONYL-SULPHIDE" ||
2981 input == "CARBONYL-SULFIDE" || input == "CARBONYL SULFIDE") {
2982 return "COS";
2983 } else if (input == "DEUT-METHANE" || input == "DEUTERIUM-METHANE" ||
2984 input == "DEUTERATED-METHANE" || input == "DEUTERATED METHANE" ||
2985 input == "DEUTERIUM METHANE" || input == "CD4") {
2986 return "CD4";
2987 } else if (input == "BF3" || input == "BORON-TRIFLUORIDE" ||
2988 input == "BORON TRIFLUORIDE") {
2989 return "BF3";
2990 } else if (input == "C2HF5" || input == "C2H2F4" || input == "C2F5H" ||
2991 input == "C2F4H2" || input == "FREON 134" || input == "FREON 134A" ||
2992 input == "FREON-134" || input == "FREON-134-A" || input == "FREON 125" ||
2993 input == "ZYRON 125" || input == "FREON-125" || input == "ZYRON-125" ||
2994 input == "TETRAFLUOROETHANE" || input == "PENTAFLUOROETHANE") {
2995 // C2H2F4 (and C2HF5).
2996 return "C2H2F4";
2997 } else if (input == "TMA" || input == "TRIMETHYLAMINE" || input == "N(CH3)3" ||
2998 input == "N-(CH3)3") {
2999 return "TMA";
3000 } else if (input == "CHF3" || input == "FREON-23" || input == "TRIFLUOROMETHANE" ||
3001 input == "FLUOROFORM") {
3002 return "CHF3";
3003 } else if (input == "CF3BR" || input == "TRIFLUOROBROMOMETHANE" ||
3004 input == "BROMOTRIFLUOROMETHANE" || input == "HALON-1301" ||
3005 input == "HALON 1301" || input == "FREON-13B1" || input == "FREON 13BI") {
3006 return "CF3Br";
3007 } else if (input == "C3F8" || input == "OCTAFLUOROPROPANE" || input == "R218" ||
3008 input == "R-218" || input == "FREON 218" || input == "FREON-218" ||
3009 input == "PERFLUOROPROPANE" || input == "RC 218" || input == "PFC 218" ||
3010 input == "RC-218" || input == "PFC-218" || input == "FLUTEC PP30" ||
3011 input == "GENETRON 218") {
3012 return "C3F8";
3013 } else if (input == "OZONE" || input == "O3") {
3014 return "O3";
3015 } else if (input == "MERCURY" || input == "HG" || input == "HG2") {
3016 return "Hg";
3017 } else if (input == "H2S" || input == "HYDROGEN SULPHIDE" || input == "SEWER GAS" ||
3018 input == "HYDROGEN-SULPHIDE" || input == "SEWER-GAS" ||
3019 input == "HYDROGEN SULFIDE" || input == "HEPATIC ACID" ||
3020 input == "HYDROGEN-SULFIDE" || input == "HEPATIC-ACID" ||
3021 input == "SULFUR HYDRIDE" || input == "DIHYDROGEN MONOSULFIDE" ||
3022 input == "SULFUR-HYDRIDE" || input == "DIHYDROGEN-MONOSULFIDE" ||
3023 input == "DIHYDROGEN MONOSULPHIDE" || input == "SULPHUR HYDRIDE" ||
3024 input == "DIHYDROGEN-MONOSULPHIDE" || input == "SULPHUR-HYDRIDE" ||
3025 input == "STINK DAMP" || input == "SULFURATED HYDROGEN" ||
3026 input == "STINK-DAMP" || input == "SULFURATED-HYDROGEN") {
3027 return "H2S";
3028 } else if (input == "N-BUTANE" || input == "N-C4H10" || input == "NBUTANE" ||
3029 input == "NC4H10") {
3030 return "nC4H10";
3031 } else if (input == "N-PENTANE" || input == "N-C5H12" || input == "NPENTANE" ||
3032 input == "NC5H12") {
3033 return "nC5H12";
3034 } else if (input == "NI-PHELPS" || input == "NI PHELPS" ||
3035 input == "NITROGEN-PHELPS" || input == "NITROGEN PHELPHS" ||
3036 input == "N2-PHELPS" || input == "N2 PHELPS" || input == "N2 (PHELPS)") {
3037 // Nitrogen
3038 return "N2 (Phelps)";
3039 } else if (input == "GERMANE" || input == "GERM" || input == "GERMANIUM-HYDRIDE" ||
3040 input == "GERMANIUM HYDRIDE" || input == "GERMANIUM TETRAHYDRIDE" ||
3041 input == "GERMANIUM-TETRAHYDRIDE" || input == "GERMANOMETHANE" ||
3042 input == "MONOGERMANE" || input == "GEH4") {
3043 return "GeH4";
3044 } else if (input == "SILANE" || input == "SIL" || input == "SILICON-HYDRIDE" ||
3045 input == "SILICON HYDRIDE" || input == "SILICON-TETRAHYDRIDE" ||
3046 input == "SILICANE" || input == "MONOSILANE" || input == "SIH4") {
3047 return "SiH4";
3048 }
3049
3050 std::cerr << m_className << "::GetGasName:\n"
3051 << " Gas " << input << " is not recognized.\n";
3052 return "";
3053}
3054
3055int MediumGas::GetGasNumberGasFile(const std::string& input) const {
3056
3057 if (input.empty()) return 0;
3058
3059 if (input == "CF4") {
3060 return 1;
3061 } else if (input == "Ar") {
3062 return 2;
3063 } else if (input == "He" || input == "He-4") {
3064 return 3;
3065 } else if (input == "He-3") {
3066 return 4;
3067 } else if (input == "Ne") {
3068 return 5;
3069 } else if (input == "Kr") {
3070 return 6;
3071 } else if (input == "Xe") {
3072 return 7;
3073 } else if (input == "CH4") {
3074 // Methane
3075 return 8;
3076 } else if (input == "C2H6") {
3077 // Ethane
3078 return 9;
3079 } else if (input == "C3H8") {
3080 // Propane
3081 return 10;
3082 } else if (input == "iC4H10") {
3083 // Isobutane
3084 return 11;
3085 } else if (input == "CO2") {
3086 return 12;
3087 } else if (input == "neoC5H12") {
3088 // Neopentane
3089 return 13;
3090 } else if (input == "H2O") {
3091 return 14;
3092 } else if (input == "O2") {
3093 return 15;
3094 } else if (input == "N2") {
3095 return 16;
3096 } else if (input == "NO") {
3097 // Nitric oxide
3098 return 17;
3099 } else if (input == "N2O") {
3100 // Nitrous oxide
3101 return 18;
3102 } else if (input == "C2H4") {
3103 // Ethene
3104 return 19;
3105 } else if (input == "C2H2") {
3106 // Acetylene
3107 return 20;
3108 } else if (input == "H2") {
3109 return 21;
3110 } else if (input == "D2") {
3111 // Deuterium
3112 return 22;
3113 } else if (input == "CO") {
3114 return 23;
3115 } else if (input == "Methylal") {
3116 // Methylal (dimethoxymethane, CH3-O-CH2-O-CH3, "hot" version)
3117 return 24;
3118 } else if (input == "DME") {
3119 return 25;
3120 } else if (input == "Reid-Step") {
3121 return 26;
3122 } else if (input == "Maxwell-Model") {
3123 return 27;
3124 } else if (input == "Reid-Ramp") {
3125 return 28;
3126 } else if (input == "C2F6") {
3127 return 29;
3128 } else if (input == "SF6") {
3129 return 30;
3130 } else if (input == "NH3") {
3131 return 31;
3132 } else if (input == "C3H6") {
3133 // Propene
3134 return 32;
3135 } else if (input == "cC3H6") {
3136 // Cyclopropane
3137 return 33;
3138 } else if (input == "CH3OH") {
3139 // Methanol
3140 return 34;
3141 } else if (input == "C2H5OH") {
3142 // Ethanol
3143 return 35;
3144 } else if (input == "C3H7OH") {
3145 // Propanol
3146 return 36;
3147 } else if (input == "Cs") {
3148 return 37;
3149 } else if (input == "F2") {
3150 // Fluorine
3151 return 38;
3152 } else if (input == "CS2") {
3153 return 39;
3154 } else if (input == "COS") {
3155 return 40;
3156 } else if (input == "CD4") {
3157 // Deuterated methane
3158 return 41;
3159 } else if (input == "BF3") {
3160 return 42;
3161 } else if (input == "C2HF5" || input == "C2H2F4") {
3162 return 43;
3163 } else if (input == "TMA") {
3164 return 44;
3165 } else if (input == "paraH2") {
3166 return 45;
3167 } else if (input == "nC3H7OH") {
3168 return 46;
3169 } else if (input == "orthoD2") {
3170 return 48;
3171 } else if (input == "CHF3") {
3172 return 50;
3173 } else if (input == "CF3Br") {
3174 return 51;
3175 } else if (input == "C3F8") {
3176 return 52;
3177 } else if (input == "O3") {
3178 // Ozone
3179 return 53;
3180 } else if (input == "Hg") {
3181 return 54;
3182 } else if (input == "H2S") {
3183 return 55;
3184 } else if (input == "nC4H10") {
3185 // n-butane
3186 return 56;
3187 } else if (input == "nC5H12") {
3188 // n-pentane
3189 return 57;
3190 } else if (input == "N2 (Phelps)") {
3191 return 58;
3192 } else if (input == "GeH4") {
3193 // Germane
3194 return 59;
3195 } else if (input == "SiH4") {
3196 // Silane
3197 return 60;
3198 }
3199
3200 std::cerr << m_className << "::GetGasNumberGasFile:\n"
3201 << " Gas " << input << " not found.\n";
3202 return 0;
3203}
3204
3205bool MediumGas::GetPhotoAbsorptionCrossSection(const double e, double& sigma,
3206 const unsigned int i) {
3207 if (i >= m_nMaxGases) {
3208 std::cerr << m_className
3209 << "::GetPhotoAbsorptionCrossSection: Index out of range.\n";
3210 return false;
3211 }
3212
3213 OpticalData optData;
3214 if (!optData.IsAvailable(m_gas[i])) return false;
3215 double eta = 0.;
3216 return optData.GetPhotoabsorptionCrossSection(m_gas[i], e, sigma, eta);
3217}
3218}
double GetMassDensity() const override
Get the mass density [g/cm3].
Definition: MediumGas.cc:299
double GetNumberDensity() const override
Get the number density [cm-3].
Definition: MediumGas.cc:293
void ReadRecord3D(std::ifstream &gasfile, double &ve, double &vb, double &vx, double &dl, double &dt, double &alpha, double &alpha0, double &eta, double &mu, double &lor, double &dis, std::array< double, 6 > &dif, std::vector< double > &rexc, std::vector< double > &rion)
Definition: MediumGas.cc:755
bool AdjustTownsendCoefficient()
Definition: MediumGas.cc:2486
std::pair< unsigned int, unsigned int > m_extrIon
Definition: MediumGas.hh:176
void GetComposition(std::string &gas1, double &f1, std::string &gas2, double &f2, std::string &gas3, double &f3, std::string &gas4, double &f4, std::string &gas5, double &f5, std::string &gas6, double &f6)
Retrieve the gas mixture.
Definition: MediumGas.cc:229
void SetAtomicNumber(const double z) override
Set the effective atomic number.
Definition: MediumGas.cc:260
bool GetGasInfo(const std::string &gasname, double &a, double &z) const
Definition: MediumGas.cc:2552
static constexpr unsigned int m_nMaxGases
Definition: MediumGas.hh:126
std::vector< std::vector< std::vector< std::vector< double > > > > m_excRates
Definition: MediumGas.hh:155
void InsertB(const int ib, const int ne, const int nb, const int na)
Definition: MediumGas.cc:1646
double m_lambdaPenningGlobal
Definition: MediumGas.hh:140
unsigned int m_intpIon
Definition: MediumGas.hh:178
std::vector< IonLevel > m_ionLevels
Definition: MediumGas.hh:172
std::array< double, m_nMaxGases > m_rPenningGas
Definition: MediumGas.hh:142
bool LoadIonMobility(const std::string &filename)
Read a table of ion mobilities as function of electric field from file.
Definition: MediumGas.cc:2215
double GetAtomicNumber() const override
Get the effective atomic number.
Definition: MediumGas.cc:303
std::array< double, m_nMaxGases > m_atNum
Definition: MediumGas.hh:132
void GetGasBits(std::bitset< 20 > &gasok) const
Definition: MediumGas.cc:2052
std::vector< std::vector< std::vector< std::vector< double > > > > m_ionRates
Definition: MediumGas.hh:156
void ZeroRowA(const int ia, const int ne, const int nb)
Definition: MediumGas.cc:1752
void ResetTables() override
Reset all tables of transport parameters.
Definition: MediumGas.cc:2292
std::vector< ExcLevel > m_excLevels
Definition: MediumGas.hh:166
void ZeroRowB(const int ib, const int ne, const int na)
Definition: MediumGas.cc:1744
MediumGas()
Constructor.
Definition: MediumGas.cc:112
int GetGasNumberGasFile(const std::string &input) const
Definition: MediumGas.cc:3055
void ReadFooter(std::ifstream &gasfile, std::array< unsigned int, 13 > &extrapH, std::array< unsigned int, 13 > &extrapL, std::array< unsigned int, 13 > &interp, unsigned int &thrAlp, unsigned int &thrAtt, unsigned int &thrDis, double &ionDiffL, double &ionDiffT, double &pgas, double &tgas)
Definition: MediumGas.cc:813
double m_temperatureTable
Definition: MediumGas.hh:149
virtual void PrintGas()
Print information about the present gas mixture and available data.
Definition: MediumGas.cc:2072
void ReadRecord1D(std::ifstream &gasfile, double &ve, double &vb, double &vx, double &dl, double &dt, double &alpha, double &alpha0, double &eta, double &mu, double &lor, double &dis, std::array< double, 6 > &dif, std::vector< double > &rexc, std::vector< double > &rion)
Definition: MediumGas.cc:789
std::array< double, m_nMaxGases > m_atWeight
Definition: MediumGas.hh:131
void GetComponent(const unsigned int i, std::string &label, double &f) override
Get the name and fraction of a given component.
Definition: MediumGas.cc:247
bool WriteGasFile(const std::string &filename)
Save the present table of gas properties (transport parameters) to a file.
Definition: MediumGas.cc:1760
bool GetMixture(const std::vector< double > &mixture, const int version, std::vector< std::string > &gasnames, std::vector< double > &percentages) const
Definition: MediumGas.cc:895
std::array< double, m_nMaxGases > m_lambdaPenningGas
Definition: MediumGas.hh:144
std::pair< unsigned int, unsigned int > m_extrExc
Definition: MediumGas.hh:175
bool GetPhotoAbsorptionCrossSection(const double e, double &sigma, const unsigned int i) override
Definition: MediumGas.cc:3205
bool MergeGasFile(const std::string &filename, const bool replaceOld)
Read table of gas properties from and merge with the existing dataset.
Definition: MediumGas.cc:932
virtual void DisablePenningTransfer()
Switch the simulation of Penning transfers off globally.
Definition: MediumGas.cc:2431
std::string GetGasName(const int gasnumber, const int version) const
Definition: MediumGas.cc:2718
double GetAtomicWeight() const override
Get the effective atomic weight.
Definition: MediumGas.cc:284
void SetNumberDensity(const double n) override
Set the number density [cm-3].
Definition: MediumGas.cc:272
bool SetComposition(const std::string &gas1, const double f1=1., const std::string &gas2="", const double f2=0., const std::string &gas3="", const double f3=0., const std::string &gas4="", const double f4=0., const std::string &gas5="", const double f5=0., const std::string &gas6="", const double f6=0.)
Set the gas mixture.
Definition: MediumGas.cc:135
void ZeroRowE(const int ie, const int nb, const int na)
Definition: MediumGas.cc:1736
void SetMassDensity(const double rho) override
Set the mass density [g/cm3].
Definition: MediumGas.cc:278
virtual bool EnablePenningTransfer(const double r, const double lambda)
Definition: MediumGas.cc:2302
std::vector< std::vector< std::vector< double > > > m_eAlp0
Definition: MediumGas.hh:152
void InsertA(const int ia, const int ne, const int nb, const int na)
Definition: MediumGas.cc:1694
bool LoadGasFile(const std::string &filename)
Read table of gas properties (transport parameters) from file.
Definition: MediumGas.cc:312
void InsertE(const int ie, const int ne, const int nb, const int na)
Definition: MediumGas.cc:1604
std::array< std::string, m_nMaxGases > m_gas
Definition: MediumGas.hh:129
bool ReadHeader(std::ifstream &gasfile, int &version, std::bitset< 20 > &gasok, bool &is3d, std::vector< double > &mixture, std::vector< double > &efields, std::vector< double > &bfields, std::vector< double > &angles, std::vector< ExcLevel > &excLevels, std::vector< IonLevel > &ionLevels)
Definition: MediumGas.cc:579
void SetAtomicWeight(const double a) override
Set the effective atomic weight.
Definition: MediumGas.cc:266
unsigned int m_intpExc
Definition: MediumGas.hh:177
std::array< double, m_nMaxGases > m_fraction
Definition: MediumGas.hh:130
Abstract base class for media.
Definition: Medium.hh:13
std::vector< double > m_bFields
Definition: Medium.hh:537
double m_pressure
Definition: Medium.hh:506
unsigned int m_intpMob
Definition: Medium.hh:591
unsigned int m_intpDis
Definition: Medium.hh:592
std::vector< std::vector< std::vector< double > > > m_eAlp
Definition: Medium.hh:546
virtual void EnableDrift(const bool on=true)
Switch electron/ion/hole on/off.
Definition: Medium.hh:67
std::string m_name
Definition: Medium.hh:502
unsigned int m_intpDif
Definition: Medium.hh:587
virtual void EnablePrimaryIonisation(const bool on=true)
Make the medium ionisable or non-ionisable.
Definition: Medium.hh:69
std::pair< unsigned int, unsigned int > m_extrAtt
Definition: Medium.hh:580
std::vector< std::vector< std::vector< double > > > m_iDifT
Definition: Medium.hh:566
std::pair< unsigned int, unsigned int > m_extrVel
Definition: Medium.hh:577
size_t SetThreshold(const std::vector< std::vector< std::vector< double > > > &tab) const
Definition: Medium.cc:1137
std::vector< std::vector< std::vector< double > > > m_eVelE
Definition: Medium.hh:541
std::vector< std::vector< std::vector< double > > > m_eVelX
Definition: Medium.hh:542
std::vector< std::vector< std::vector< double > > > m_eDifL
Definition: Medium.hh:544
void Init(const size_t nE, const size_t nB, const size_t nA, std::vector< std::vector< std::vector< double > > > &tab, const double val)
Definition: Medium.cc:1295
unsigned int m_intpVel
Definition: Medium.hh:586
std::vector< double > m_eFields
Definition: Medium.hh:536
std::vector< std::vector< std::vector< double > > > m_eAtt
Definition: Medium.hh:547
std::vector< double > m_bAngles
Definition: Medium.hh:538
std::vector< std::vector< std::vector< double > > > m_iDifL
Definition: Medium.hh:565
std::vector< std::vector< std::vector< double > > > m_eLor
Definition: Medium.hh:548
std::vector< std::vector< std::vector< double > > > m_eDifT
Definition: Medium.hh:545
unsigned int m_intpAtt
Definition: Medium.hh:589
std::pair< unsigned int, unsigned int > m_extrDis
Definition: Medium.hh:583
std::pair< unsigned int, unsigned int > m_extrAlp
Definition: Medium.hh:579
std::vector< std::vector< std::vector< double > > > m_eVelB
Definition: Medium.hh:543
unsigned int m_eThrAtt
Definition: Medium.hh:571
std::pair< unsigned int, unsigned int > m_extrLor
Definition: Medium.hh:581
unsigned int m_nComponents
Definition: Medium.hh:500
unsigned int m_intpLor
Definition: Medium.hh:590
std::string m_className
Definition: Medium.hh:493
std::pair< unsigned int, unsigned int > m_extrDif
Definition: Medium.hh:578
unsigned int m_iThrDis
Definition: Medium.hh:574
virtual void ResetTables()
Reset all tables of transport parameters.
Definition: Medium.cc:901
std::pair< unsigned int, unsigned int > m_extrMob
Definition: Medium.hh:582
std::vector< std::vector< std::vector< std::vector< double > > > > m_eDifM
Definition: Medium.hh:550
bool m_isChanged
Definition: Medium.hh:527
std::vector< std::vector< std::vector< double > > > m_iMob
Definition: Medium.hh:564
unsigned int m_eThrAlp
Definition: Medium.hh:570
double m_temperature
Definition: Medium.hh:504
std::vector< std::vector< std::vector< double > > > m_iDis
Definition: Medium.hh:567
bool SetIonMobility(const std::vector< double > &fields, const std::vector< double > &mobilities)
Definition: Medium.cc:1039
unsigned int m_intpAlp
Definition: Medium.hh:588
Photoabsorption cross-sections for some gases.
Definition: OpticalData.hh:11
bool GetPhotoabsorptionCrossSection(const std::string &material, const double e, double &cs, double &eta)
Photo-absorption cross-section and ionisation yield at a given energy.
Definition: OpticalData.cc:46
bool IsAvailable(const std::string &material) const
Check whether optical data have been implemented for a given gas.
Definition: OpticalData.cc:39
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615