31void ClearBank(std::vector<Heed::gparticle*>& bank) {
32 for (
auto particle : bank)
33 if (particle)
delete particle;
51 m_conductionElectrons.reserve(1000);
52 m_conductionIons.reserve(1000);
60 const double t0,
const double dx0,
const double dy0,
62 m_hasActiveTrack =
false;
67 std::cerr <<
m_className <<
"::NewTrack: Sensor is not defined.\n";
72 if (!UpdateBoundingBox(update))
return false;
78 <<
" No medium at initial position.\n";
82 <<
" Medium at initial position is not ionisable.\n";
87 if (medium->
GetName() != m_mediumName ||
95 if (!Setup(medium))
return false;
97 m_mediumName = medium->
GetName();
103 m_deltaElectrons.clear();
104 m_conductionElectrons.clear();
105 m_conductionIons.clear();
108 double dx = dx0, dy = dy0, dz = dz0;
109 const double d = sqrt(dx * dx + dy * dy + dz * dz);
113 <<
" Direction vector has zero norm.\n"
114 <<
" Initial direction is randomized.\n";
125 velocity = velocity * Heed::CLHEP::c_light *
GetBeta();
128 std::cout <<
m_className <<
"::NewTrack:\n Track starts at (" << x0
129 <<
", " << y0 <<
", " << z0 <<
") at time " << t0 <<
"\n"
130 <<
" Direction: (" << dx <<
", " << dy <<
", " << dz <<
")\n";
135 Heed::point p0((x0 - m_cX) * 10., (y0 - m_cY) * 10., (z0 - m_cZ) * 10.);
181 m_radStraight * Heed::CLHEP::cm,
182 m_stepAngleStraight * Heed::CLHEP::rad,
183 m_stepAngleCurved * Heed::CLHEP::rad);
185 particle.
fly(m_particleBank);
186 m_bankIterator = m_particleBank.begin();
187 m_hasActiveTrack =
true;
197 std::cerr <<
m_className <<
"::GetClusterDensity:\n"
198 <<
" Track has not been initialized.\n";
203 std::cerr <<
m_className <<
"::GetClusterDensity:\n"
204 <<
" Ionisation cross-section is not available.\n";
208 return m_transferCs->quanC;
213 std::cerr <<
m_className <<
"::GetStoppingPower:\n"
214 <<
" Track has not been initialized.\n";
219 std::cerr <<
m_className <<
"::GetStoppingPower:\n"
220 <<
" Ionisation cross-section is not available.\n";
224 return m_transferCs->meanC1 * 1.e6;
228 double& tcls,
int& n,
double& e,
double& extra) {
230 return GetCluster(xcls, ycls, zcls, tcls, n, ni, e, extra);
234 double& tcls,
int& ne,
int& ni,
double& e,
237 xcls = ycls = zcls = tcls = 0.;
242 m_deltaElectrons.clear();
243 m_conductionElectrons.clear();
244 m_conductionIons.clear();
249 <<
" Track has not been initialized. Call NewTrack first.\n";
253 if (m_particleBank.empty())
return false;
254 std::vector<Heed::gparticle*>::const_iterator end = m_particleBank.end();
255 if (m_bankIterator == end)
return false;
259 for (; m_bankIterator != end; ++m_bankIterator) {
262 if (!virtualPhoton) {
264 <<
" Particle is not a virtual photon. Program bug!\n";
271 xcls = virtualPhoton->
position().
x * 0.1 + m_cX;
272 ycls = virtualPhoton->
position().
y * 0.1 + m_cY;
273 zcls = virtualPhoton->
position().
z * 0.1 + m_cZ;
274 tcls = virtualPhoton->
time();
276 if (!IsInside(xcls, ycls, zcls))
continue;
278 m_conductionIons.emplace_back(
285 if (!virtualPhoton)
return false;
289 std::vector<Heed::gparticle*> secondaries;
291 virtualPhoton->
fly(secondaries);
295 while (!secondaries.empty()) {
296 std::vector<Heed::gparticle*> newSecondaries;
298 for (
auto secondary : secondaries) {
303 const double x = delta->position().x * 0.1 + m_cX;
304 const double y = delta->position().y * 0.1 + m_cY;
305 const double z = delta->position().z * 0.1 + m_cZ;
306 if (!IsInside(x, y, z))
continue;
307 if (m_doDeltaTransport) {
309 delta->fly(newSecondaries);
311 m_conductionElectrons.insert(m_conductionElectrons.end(),
312 delta->conduction_electrons.begin(),
313 delta->conduction_electrons.end());
314 m_conductionIons.insert(m_conductionIons.end(),
315 delta->conduction_ions.begin(),
316 delta->conduction_ions.end());
319 deltaElectron newDeltaElectron;
320 newDeltaElectron.x = delta->position().x * 0.1 + m_cX;
321 newDeltaElectron.y = delta->position().y * 0.1 + m_cY;
322 newDeltaElectron.z = delta->position().z * 0.1 + m_cZ;
323 newDeltaElectron.t = delta->time();
324 newDeltaElectron.e = delta->kinetic_energy() * 1.e6;
325 newDeltaElectron.dx = delta->direction().x;
326 newDeltaElectron.dy = delta->direction().y;
327 newDeltaElectron.dz = delta->direction().z;
328 m_deltaElectrons.push_back(std::move(newDeltaElectron));
336 <<
" Particle is neither an electron nor a photon.\n";
338 extra += photon->m_energy * 1.e6;
339 const double x = photon->position().x * 0.1 + m_cX;
340 const double y = photon->position().y * 0.1 + m_cY;
341 const double z = photon->position().z * 0.1 + m_cZ;
342 if (!IsInside(x, y, z))
continue;
344 if (m_usePhotonReabsorption) photon->fly(newSecondaries);
346 for (
auto secondary : secondaries)
347 if (secondary)
delete secondary;
349 secondaries.swap(newSecondaries);
352 ne = m_doDeltaTransport ? m_conductionElectrons.size()
353 : m_deltaElectrons.size();
354 ni = m_conductionIons.size();
359 double& z,
double& t,
double& e,
double& dx,
360 double& dy,
double& dz) {
364 <<
" Track has not been initialized. Call NewTrack first.\n";
368 if (m_doDeltaTransport) {
370 if (i >= m_conductionElectrons.size()) {
371 std::cerr <<
m_className <<
"::GetElectron: Index out of range.\n";
375 x = m_conductionElectrons[i].x * 0.1 + m_cX;
376 y = m_conductionElectrons[i].y * 0.1 + m_cY;
377 z = m_conductionElectrons[i].z * 0.1 + m_cZ;
378 t = m_conductionElectrons[i].time;
384 if (i >= m_deltaElectrons.size()) {
386 <<
" Delta electron number out of range.\n";
390 x = m_deltaElectrons[i].x;
391 y = m_deltaElectrons[i].y;
392 z = m_deltaElectrons[i].z;
393 t = m_deltaElectrons[i].t;
394 e = m_deltaElectrons[i].e;
395 dx = m_deltaElectrons[i].dx;
396 dy = m_deltaElectrons[i].dy;
397 dz = m_deltaElectrons[i].dz;
405 if (i >= m_conductionIons.size()) {
406 std::cerr <<
m_className <<
"::GetIon: Index out of range.\n";
410 x = m_conductionIons[i].x * 0.1 + m_cX;
411 y = m_conductionIons[i].y * 0.1 + m_cY;
412 z = m_conductionIons[i].z * 0.1 + m_cZ;
413 t = m_conductionIons[i].time;
418 const double z0,
const double t0,
419 const double e0,
const double dx0,
420 const double dy0,
const double dz0,
427 const double z0,
const double t0,
428 const double e0,
const double dx0,
429 const double dy0,
const double dz0,
435 if (!m_doDeltaTransport) {
436 std::cerr <<
m_className <<
"::TransportDeltaElectron:\n"
437 <<
" Delta electron transport has been switched off.\n";
443 std::cerr <<
m_className <<
"::TransportDeltaElectron:\n"
444 <<
" Sensor is not defined.\n";
450 if (!UpdateBoundingBox(update))
return;
455 std::cerr <<
m_className <<
"::TransportDeltaElectron:\n"
456 <<
" No medium at initial position.\n";
459 std::cerr <<
"TrackHeed:TransportDeltaElectron:\n"
460 <<
" Medium at initial position is not ionisable.\n";
466 if (medium->
GetName() != m_mediumName ||
471 m_hasActiveTrack =
false;
476 if (!Setup(medium))
return;
478 m_mediumName = medium->
GetName();
482 m_deltaElectrons.clear();
483 m_conductionElectrons.clear();
484 m_conductionIons.clear();
488 Heed::point p0((x0 - m_cX) * 10., (y0 - m_cY) * 10., (z0 - m_cZ) * 10.);
499 double dx = dx0, dy = dy0, dz = dz0;
500 const double d = sqrt(dx * dx + dy * dy + dz * dz);
513 const double gamma = 1. + e0 / ElectronMass;
514 const double beta = sqrt(1. - 1. / (gamma * gamma));
515 double speed = Heed::CLHEP::c_light * beta;
516 velocity = velocity * speed;
519 std::vector<Heed::gparticle*> secondaries;
522 delta.
fly(secondaries);
523 ClearBank(secondaries);
527 nel = m_conductionElectrons.size();
528 ni = m_conductionIons.size();
532 const double z0,
const double t0,
533 const double e0,
const double dx0,
534 const double dy0,
const double dz0,
int& nel) {
540 const double z0,
const double t0,
541 const double e0,
const double dx0,
542 const double dy0,
const double dz0,
int& nel,
550 <<
" Photon energy must be positive.\n";
556 std::cerr <<
m_className <<
"::TransportPhoton: Sensor is not defined.\n";
562 if (!UpdateBoundingBox(update))
return;
568 <<
" No medium at initial position.\n";
571 std::cerr <<
"TrackHeed:TransportPhoton:\n"
572 <<
" Medium at initial position is not ionisable.\n";
578 if (medium->
GetName() != m_mediumName ||
587 if (!Setup(medium))
return;
589 m_mediumName = medium->
GetName();
595 m_hasActiveTrack =
false;
597 m_deltaElectrons.clear();
598 m_conductionElectrons.clear();
599 m_conductionIons.clear();
602 double dx = dx0, dy = dy0, dz = dz0;
603 const double d = sqrt(dx * dx + dy * dy + dz * dz);
614 velocity = velocity * Heed::CLHEP::c_light;
618 Heed::point p0((x0 - m_cX) * 10., (y0 - m_cY) * 10., (z0 - m_cZ) * 10.);
623 std::vector<Heed::gparticle*> secondaries;
624 photon.
fly(secondaries);
626 while (!secondaries.empty()) {
627 std::vector<Heed::gparticle*> newSecondaries;
629 std::vector<Heed::gparticle*>::iterator it;
630 for (it = secondaries.begin(); it != secondaries.end(); ++it) {
634 if (m_doDeltaTransport) {
636 delta->
fly(newSecondaries);
638 m_conductionElectrons.insert(m_conductionElectrons.end(),
639 delta->conduction_electrons.begin(),
640 delta->conduction_electrons.end());
641 m_conductionIons.insert(m_conductionIons.end(),
642 delta->conduction_ions.begin(),
643 delta->conduction_ions.end());
646 deltaElectron newDeltaElectron;
647 newDeltaElectron.x = delta->position().x * 0.1 + m_cX;
648 newDeltaElectron.y = delta->position().y * 0.1 + m_cY;
649 newDeltaElectron.z = delta->position().z * 0.1 + m_cZ;
650 newDeltaElectron.t = delta->time();
651 newDeltaElectron.e = delta->kinetic_energy() * 1.e6;
652 newDeltaElectron.dx = delta->direction().x;
653 newDeltaElectron.dy = delta->direction().y;
654 newDeltaElectron.dz = delta->direction().z;
655 m_deltaElectrons.push_back(std::move(newDeltaElectron));
661 if (!fluorescencePhoton) {
663 <<
" Unknown secondary particle.\n";
664 ClearBank(secondaries);
665 ClearBank(newSecondaries);
668 fluorescencePhoton->fly(newSecondaries);
670 secondaries.swap(newSecondaries);
671 ClearBank(newSecondaries);
673 ClearBank(secondaries);
675 nel = m_doDeltaTransport ? m_conductionElectrons.size()
676 : m_deltaElectrons.size();
677 ni = m_conductionIons.size();
687 if (fabs(e1 - e0) < Small) {
689 <<
" Invalid energy range:\n"
690 <<
" " << e0 <<
" < E [eV] < " << e1 <<
"\n";
696 <<
" Number of intervals must be > 0.\n";
700 m_emin = std::min(e0, e1);
701 m_emax = std::max(e0, e1);
704 m_nEnergyIntervals = nsteps;
708 if (fabs(z) < Small) {
710 <<
" Particle cannot have zero charge.\n";
715 <<
" Particle mass must be greater than zero.\n";
724bool TrackHeed::Setup(
Medium* medium) {
726 std::string databasePath;
727 char* dbPath = std::getenv(
"HEED_DATABASE");
728 if (dbPath == NULL) {
730 dbPath = std::getenv(
"GARFIELD_HOME");
731 if (dbPath == NULL) {
732 std::cerr <<
m_className <<
"::Setup:\n Cannot retrieve database path "
733 <<
"(environment variables HEED_DATABASE and GARFIELD_HOME "
734 <<
"are not defined).\n Cannot proceed.\n";
737 databasePath = std::string(dbPath) +
"/Heed/heed++/database";
739 databasePath = dbPath;
741 if (databasePath[databasePath.size() - 1] !=
'/') {
742 databasePath.append(
"/");
747 std::cerr <<
m_className <<
"::Setup: Null pointer.\n";
752 m_energyMesh.reset(
new Heed::EnergyMesh(m_emin, m_emax, m_nEnergyIntervals));
754 if (medium->
IsGas()) {
755 if (!SetupGas(medium))
return false;
757 if (!SetupMaterial(medium))
return false;
766 if (!SetupDelta(databasePath))
return false;
769 const double nc = m_transferCs->quanC;
770 const double dedx = m_transferCs->meanC * 1.e3;
771 const double dedx1 = m_transferCs->meanC1 * 1.e3;
772 const double w = m_matter->W * 1.e6;
773 const double f = m_matter->F;
774 const double minI = m_matter->min_ioniz_pot * 1.e6;
776 std::cout <<
" Cluster density: " << nc <<
" cm-1\n";
777 std::cout <<
" Stopping power (restricted): " << dedx <<
" keV/cm\n";
778 std::cout <<
" Stopping power (incl. tail): " << dedx1 <<
" keV/cm\n";
779 std::cout <<
" W value: " << w <<
" eV\n";
780 std::cout <<
" Fano factor: " << f <<
"\n";
781 std::cout <<
" Min. ionization potential: " << minI <<
" eV\n";
786 m_chamber.reset(
new HeedChamber(primSys, m_lX, m_lY, m_lZ,
787 *m_transferCs.get(), *m_deltaCs.get()));
792bool TrackHeed::SetupGas(Medium* medium) {
794 double pressure = medium->GetPressure();
795 pressure = (pressure / AtmosphericPressure) * Heed::CLHEP::atmosphere;
796 double temperature = medium->GetTemperature();
798 const int nComponents = medium->GetNumberOfComponents();
799 if (nComponents < 1) {
801 std::cerr <<
" Gas " << medium->GetName() <<
" has zero constituents.\n";
805 std::vector<Heed::MolecPhotoAbsCS*> molPacs(nComponents,
nullptr);
806 std::vector<std::string> notations;
807 std::vector<double> fractions;
809 for (
int i = 0; i < nComponents; ++i) {
812 medium->GetComponent(i, gasname, frac);
814 if (gasname ==
"He-3") gasname =
"He";
815 if (gasname ==
"CD4") gasname =
"CH4";
816 if (gasname ==
"iC4H10" || gasname ==
"nC4H10") gasname =
"C4H10";
817 if (gasname ==
"neoC5H12" || gasname ==
"nC5H12") gasname =
"C5H12";
818 if (gasname ==
"H2O") gasname =
"Water";
819 if (gasname ==
"D2") gasname =
"H2";
820 if (gasname ==
"cC3H6") gasname =
"C3H6";
822 if (gasname ==
"CF4")
824 else if (gasname ==
"Ar")
826 else if (gasname ==
"He")
828 else if (gasname ==
"Ne")
830 else if (gasname ==
"Kr")
832 else if (gasname ==
"Xe")
834 else if (gasname ==
"CH4")
836 else if (gasname ==
"C2H6")
838 else if (gasname ==
"C3H8")
840 else if (gasname ==
"C4H10")
842 else if (gasname ==
"CO2")
844 else if (gasname ==
"C5H12")
846 else if (gasname ==
"Water")
848 else if (gasname ==
"O2")
850 else if (gasname ==
"N2" || gasname ==
"N2 (Phelps)")
852 else if (gasname ==
"NO")
854 else if (gasname ==
"N2O")
856 else if (gasname ==
"C2H4")
858 else if (gasname ==
"C2H2")
860 else if (gasname ==
"H2")
862 else if (gasname ==
"CO")
864 else if (gasname ==
"Methylal")
866 else if (gasname ==
"DME")
868 else if (gasname ==
"C2F6")
870 else if (gasname ==
"SF6")
872 else if (gasname ==
"NH3")
874 else if (gasname ==
"C3H6")
876 else if (gasname ==
"CH3OH")
878 else if (gasname ==
"C2H5OH")
880 else if (gasname ==
"C3H7OH")
882 else if (gasname ==
"Cs")
884 else if (gasname ==
"F2")
886 else if (gasname ==
"CS2")
888 else if (gasname ==
"COS")
890 else if (gasname ==
"CD4")
892 else if (gasname ==
"BF3")
894 else if (gasname ==
"C2HF5")
896 else if (gasname ==
"C2H2F4")
898 else if (gasname ==
"CHF3")
900 else if (gasname ==
"CF3Br")
902 else if (gasname ==
"C3F8")
904 else if (gasname ==
"O3")
906 else if (gasname ==
"Hg")
908 else if (gasname ==
"H2S")
910 else if (gasname ==
"GeH4")
912 else if (gasname ==
"SiH4")
915 std::cerr <<
m_className <<
"::SetupGas:\n Photoabsorption "
916 <<
"cross-section for " << gasname <<
" not available.\n";
919 notations.push_back(gasname);
920 fractions.push_back(frac);
922 if (m_usePacsOutput) {
923 std::ofstream pacsfile;
924 pacsfile.open(
"heed_pacs.txt", std::ios::out);
925 const int nValues = m_energyMesh->get_q();
927 for (
int i = 0; i < nValues; ++i) {
928 double e = m_energyMesh->get_e(i);
929 pacsfile << 1.e6 * e <<
" ";
930 for (
int j = 0; j < nComponents; ++j) {
931 pacsfile << molPacs[j]->get_ACS(e) <<
" " << molPacs[j]->get_ICS(e)
940 const std::string gasname = FindUnusedMaterialName(medium->GetName());
941 m_gas.reset(
new Heed::GasDef(gasname, gasname, nComponents, notations,
942 fractions, pressure, temperature, -1.));
944 const double w = std::max(medium->GetW() * 1.e-6, 0.);
945 double f = medium->GetFanoFactor();
954bool TrackHeed::SetupMaterial(Medium* medium) {
956 double temperature = medium->GetTemperature();
957 const double density =
958 medium->GetMassDensity() * Heed::CLHEP::gram / Heed::CLHEP::cm3;
960 const int nComponents = medium->GetNumberOfComponents();
961 std::vector<Heed::AtomPhotoAbsCS*> atPacs(nComponents,
nullptr);
963 std::vector<std::string> notations;
964 std::vector<double> fractions;
965 for (
int i = 0; i < nComponents; ++i) {
966 std::string materialName;
968 medium->GetComponent(i, materialName, frac);
969 if (materialName ==
"C")
971 else if (materialName ==
"Si")
974 else if (materialName ==
"Ga")
976 else if (materialName ==
"Ge")
978 else if (materialName ==
"As")
980 else if (materialName ==
"Cd")
982 else if (materialName ==
"Te")
986 std::cerr <<
" Photoabsorption cross-section data for " << materialName
987 <<
" are not implemented.\n";
990 notations.push_back(materialName);
991 fractions.push_back(frac);
993 if (m_usePacsOutput) {
994 std::ofstream pacsfile;
995 pacsfile.open(
"heed_pacs.txt", std::ios::out);
996 const int nValues = m_energyMesh->get_q();
998 for (
int i = 0; i < nValues; ++i) {
999 double e = m_energyMesh->get_e(i);
1000 pacsfile << 1.e6 * e <<
" ";
1001 for (
int j = 0; j < nComponents; ++j) {
1002 pacsfile << atPacs[j]->get_ACS(e) <<
" " << atPacs[j]->get_ICS(e)
1010 const std::string materialName = FindUnusedMaterialName(medium->GetName());
1011 m_material.reset(
new Heed::MatterDef(materialName, materialName, nComponents,
1012 notations, fractions, density,
1015 double w = medium->GetW() * 1.e-6;
1017 double f = medium->GetFanoFactor();
1026bool TrackHeed::SetupDelta(
const std::string& databasePath) {
1028 std::string filename = databasePath +
"cbdel.dat";
1031 filename = databasePath +
"elastic_disp.dat";
1036 const double w = m_matter->W * 1.e6;
1037 const double f = m_matter->F;
1038 filename = databasePath +
"delta_path.dat";
1042 m_matter.get(), m_elScat.get(), m_lowSigma.get(), m_pairProd.get()));
1049std::string TrackHeed::FindUnusedMaterialName(
const std::string& namein) {
1050 std::string nameout = namein;
1051 unsigned int counter = 0;
1053 nameout = namein +
"_" + std::to_string(counter);
1059void TrackHeed::ClearParticleBank() {
1061 ClearBank(m_particleBank);
1062 m_bankIterator = m_particleBank.end();
1065bool TrackHeed::IsInside(
const double x,
const double y,
const double z) {
1069 Medium* medium =
nullptr;
1072 if (medium->GetName() != m_mediumName ||
1073 fabs(medium->GetMassDensity() - m_mediumDensity) > 1.e-9 ||
1074 !medium->IsIonisable()) {
1080bool TrackHeed::UpdateBoundingBox(
bool& update) {
1082 double xmin = 0., ymin = 0., zmin = 0.;
1083 double xmax = 0., ymax = 0., zmax = 0.;
1085 std::cerr <<
m_className <<
"::UpdateBoundingBox: Drift area is not set.\n";
1090 const double lx =
fabs(xmax - xmin);
1091 const double ly =
fabs(ymax - ymin);
1092 const double lz =
fabs(zmax - zmin);
1094 std::cout <<
m_className <<
"::UpdateBoundingBox:\n"
1095 <<
" Bounding box dimensions:\n"
1096 <<
" x: " << lx <<
" cm\n"
1097 <<
" y: " << ly <<
" cm\n"
1098 <<
" z: " << lz <<
" cm\n";
1100 if (
fabs(lx - m_lX) > Small ||
fabs(ly - m_lY) > Small ||
1101 fabs(lz - m_lZ) > Small) {
1107 m_hasActiveTrack =
false;
1110 m_cX = (std::isinf(xmin) || std::isinf(xmax)) ? 0. : 0.5 * (xmin + xmax);
1111 m_cY = (std::isinf(ymin) || std::isinf(ymax)) ? 0. : 0.5 * (ymin + ymax);
1112 m_cZ = (std::isinf(zmin) || std::isinf(zmax)) ? 0. : 0.5 * (zmin + zmax);
1114 std::cout <<
m_className <<
"::UpdateBoundingBox:\n"
1115 <<
" Center of bounding box:\n"
1116 <<
" x: " << m_cX <<
" cm\n"
1117 <<
" y: " << m_cY <<
" cm\n"
1118 <<
" z: " << m_cZ <<
" cm\n";
1122 m_fieldMap->SetCentre(m_cX, m_cY, m_cZ);
Abstract base class for media.
virtual double GetMassDensity() const
Get the mass density [g/cm3].
const std::string & GetName() const
Get the medium name/identifier.
virtual bool IsGas() const
Is this medium a gas?
bool IsIonisable() const
Is charge deposition by charged particles/photon enabled in this medium?
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
bool GetArea(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Return the current user area.
void EnableMagneticField()
Take the magnetic field into account in the stepping algorithm.
void SetEnergyMesh(const double e0, const double e1, const int nsteps)
bool GetCluster(double &xcls, double &ycls, double &zcls, double &tcls, int &n, double &e, double &extra) override
void EnableElectricField()
Take the electric field into account in the stepping algorithm.
void DisableMagneticField()
Do not take the magnetic field into account in the stepping algorithm.
void TransportDeltaElectron(const double x0, const double y0, const double z0, const double t0, const double e0, const double dx0, const double dy0, const double dz0, int &ne, int &ni)
double GetFanoFactor() const
Return the Fano factor of the medium (of the last simulated track).
bool NewTrack(const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0) override
bool GetElectron(const unsigned int i, double &x, double &y, double &z, double &t, double &e, double &dx, double &dy, double &dz)
void TransportPhoton(const double x0, const double y0, const double z0, const double t0, const double e0, const double dx0, const double dy0, const double dz0, int &ne, int &ni)
bool GetIon(const unsigned int i, double &x, double &y, double &z, double &t) const
virtual ~TrackHeed()
Destructor.
double GetW() const
Return the W value of the medium (of the last simulated track).
void SetParticleUser(const double m, const double z)
double GetClusterDensity() override
void DisableElectricField()
Do not take the electric field into account in the stepping algorithm.
double GetStoppingPower() override
Get the stopping power (mean energy loss [eV] per cm).
Abstract base class for track generation.
double GetBeta() const
Return the speed ( ) of the projectile.
double GetGamma() const
Return the Lorentz factor of the projectile.
std::string m_particleName
void PlotCluster(const double x0, const double y0, const double z0)
void PlotNewTrack(const double x0, const double y0, const double z0)
std::vector< HeedCondElectron > conduction_electrons
std::vector< HeedCondElectron > conduction_ions
Retrieve electric and magnetic field from Sensor.
double m_energy
Photon energy [MeV].
static MatterDef * get_MatterDef(const std::string &fnotation)
vfloat time() const
Get the current time of the particle.
const vec & position() const
Get the current position of the particle.
virtual void fly(std::vector< gparticle * > &secondaries)
Transport the particle.
void set_step_limits(const vfloat fmax_range, const vfloat frad_for_straight, const vfloat fmax_straight_arange, const vfloat fmax_circ_arange)
Set limits/parameters for trajectory steps.
double kinetic_energy() const
Get the current kinetic energy.
void set_mass(const double m)
void set_charge(const double z)
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
ExAtomPhotoAbsCS Gallium_for_GaAs_PACS(31, shelllist_dir_name+"shelllist_solid.dat", pacs_table_dir_name+"Ga.dat", "Ga_for_GaAs")
ExAtomPhotoAbsCS Tellurium_for_CdTe_PACS(52, shelllist_dir_name+"shelllist_solid.dat", pacs_table_dir_name+"Te.dat", "Te_for_CdTe")
MolecPhotoAbsCS C2H2_MPACS
MolecPhotoAbsCS Methylal_MPACS
particle_def pi_minus_meson_def("pi_minus_meson", "pi-", 139.56755 *MeV/c_squared, -eplus, 0, 0, 0.0, spin_def(1.0, -1.0))
MolecPhotoAbsCS CH3OH_MPACS
particle_def pi_plus_meson_def("pi_plus_meson", "pi+", 139.56755 *MeV/c_squared, eplus, 0, 0, 0.0, spin_def(1.0, 1.0))
MolecPhotoAbsCS CF4_MPACS
ExAtomPhotoAbsCS Silicon_crystal_PACS(14, shelllist_dir_name+"shelllist_solid.dat", pacs_table_dir_name+"Si.dat", "Si_crystal")
particle_def muon_minus_def("muon_minus", "mu-", 105.658367 *MeV/c_squared, electron_charge, 1, 0, 0.5, spin_def(0.0, 0.0))
MolecPhotoAbsCS CS2_MPACS
MolecPhotoAbsCS NH3_MPACS
MolecPhotoAbsCS CF3Br_MPACS
long last_particle_number
particle_def alpha_particle_def("alpha_particle", "alpha", 3727.417 *MeV/c_squared, 2 *eplus, 0, 4, 0.0, spin_def(0.0, 0.0))
particle_def anti_proton_def("", "p-", proton_def)
MolecPhotoAbsCS C2H6_MPACS
particle_def deuteron_def("deuteron", "dtr", 1875.613 *MeV/c_squared, eplus, 0, 2, 0.0, spin_def(0.0, 0.0))
MolecPhotoAbsCS SiH4_MPACS
particle_def proton_def("proton", "p+", proton_mass_c2/c_squared, eplus, 0, 1, 0.5, spin_def(0.5, 0.5))
MolecPhotoAbsCS SF6_MPACS
MolecPhotoAbsCS C2H2F4_MPACS
MolecPhotoAbsCS H2O_MPACS
MolecPhotoAbsCS C2H4_MPACS
MolecPhotoAbsCS C2F6_MPACS
MolecPhotoAbsCS CO2_MPACS
MolecPhotoAbsCS COS_MPACS
particle_def K_minus_meson_def("K_minus_meson_def", "K-", K_plus_meson_def)
MolecPhotoAbsCS N2O_MPACS
ExAtomPhotoAbsCS Arsenic_for_GaAs_PACS(33, shelllist_dir_name+"shelllist_solid.dat", pacs_table_dir_name+"As.dat", "As_for_GaAs")
particle_def K_plus_meson_def("K_plus_meson_def", "K+", 493.677 *MeV/c_squared, 1, 0, 0, 0.0, spin_def(0.5, -0.5))
MolecPhotoAbsCS C2HF5_MPACS
particle_def user_particle_def("user_particle", "X", 139.56755 *MeV/c_squared, eplus, 0, 0, 0.0, spin_def(0.0, 0.0))
ExAtomPhotoAbsCS Germanium_PACS(32, shelllist_dir_name+"shelllist.dat", pacs_table_dir_name+"Ge.dat")
ExAtomPhotoAbsCS Cadmium_for_CdTe_PACS(48, shelllist_dir_name+"shelllist_solid.dat", pacs_table_dir_name+"Cd.dat", "Cd_for_CdTe")
DoubleAc fabs(const DoubleAc &f)
MolecPhotoAbsCS BF3_MPACS
MolecPhotoAbsCS H2S_MPACS
particle_def muon_plus_def("muon_plus", "mu+", muon_minus_def)
particle_def positron_def("positron", "e+", electron_def)
MolecPhotoAbsCS DME_MPACS
particle_def electron_def("electron", "e-", electron_mass_c2/c_squared, electron_charge, 1, 0, 0.5, spin_def(0.0, 0.0))
MolecPhotoAbsCS C3H7OH_MPACS
MolecPhotoAbsCS C3F8_MPACS
MolecPhotoAbsCS C3H6_MPACS
MolecPhotoAbsCS C4H10_MPACS
ExAtomPhotoAbsCS Carbon_PACS(6, shelllist_dir_name+"shelllist.dat", pacs_table_dir_name+"C.dat")
MolecPhotoAbsCS C3H8_MPACS
constexpr double standard_factor_Fano
MolecPhotoAbsCS C5H12_MPACS
MolecPhotoAbsCS C2H5OH_MPACS
MolecPhotoAbsCS CHF3_MPACS
MolecPhotoAbsCS CH4_MPACS
MolecPhotoAbsCS GeH4_MPACS