MIDSX 0.1
A x-ray transport code system for dosimetry
Loading...
Searching...
No Matches
material_data.h
1#ifndef MCXRAYTRANSPORT_MATERIAL_DATA_H
2#define MCXRAYTRANSPORT_MATERIAL_DATA_H
3
4#include "data_access_object.h"
5#include "material_properties.h"
6#include "material_helpers.h"
7#include "interpolators.h"
8#include "probability_dist.h"
9#include "constants.h"
10#include <string>
11#include <memory>
12#include <Eigen/Dense>
13
20public:
28
37 Eigen::Matrix<double, Eigen::Dynamic, 2> getIncoherentScatteringCrossSectionMatrix() const { return incoherent_cs_matrix_; }
38
47 Eigen::Matrix<double, Eigen::Dynamic, 2> getCoherentScatteringCrossSectionMatrix() const { return coherent_cs_matrix_; }
48
57 Eigen::Matrix<double, Eigen::Dynamic, 2> getPhotoelectricCrossSectionMatrix() const { return photoelectric_cs_matrix_; }
58
67 Eigen::Matrix<double, Eigen::Dynamic, 2> getTotalCrossSectionMatrix() const { return total_cs_matrix_; }
68
77 Eigen::Matrix<double, Eigen::Dynamic, 2> getIncoherentScatteringFunctionMatrix() const { return incoherent_scattering_function_matrix_; }
78
87 Eigen::Matrix<double, Eigen::Dynamic, 2> getCoherentFormFactorMatrix() const { return coherent_form_factor_matrix_; }
88
95 double interpolateIncoherentScatteringCrossSection(double energy) const { return (incoherent_cs_interpolator_)(energy); }
96
103 double interpolateCoherentScatteringCrossSection(double energy) const { return (coherent_cs_interpolator_)(energy); }
104
111 double interpolatePhotoelectricCrossSection(double energy) const { return (photoelectric_cs_interpolator_)(energy); }
112
119 double interpolateTotalCrossSection(double energy) const { return (total_cs_interpolator_)(energy); }
120
127 double interpolateIncoherentScatteringFunction(double x) const { return (incoherent_scattering_function_interpolator_)(x); }
128
135 double interpolateCoherentFormFactor(double x) const { return (coherent_form_factor_interpolator_)(x); }
136
143 double interpolateMassEnergyAbsorptionCoefficient(double energy) const { return (mass_energy_absorption_coefficient_interpolator_)(energy); }
144
151 double sampleCoherentScatteringDCS(double energy) const { return coherent_scattering_dcs_dist_.sample(energy); }
152
153private:
154 MaterialProperties& properties_;
155 DataAccessObject& dao_;
156
157 Eigen::Matrix<double, Eigen::Dynamic, 2> incoherent_cs_matrix_;
158 Interpolator::LogLogSpline incoherent_cs_interpolator_;
159 Eigen::Matrix<double, Eigen::Dynamic, 2> coherent_cs_matrix_;
160 Interpolator::LogLogSpline coherent_cs_interpolator_;
161 Eigen::Matrix<double, Eigen::Dynamic, 2> photoelectric_cs_matrix_;
162 Interpolator::LogLogLinear photoelectric_cs_interpolator_;
163 Eigen::Matrix<double, Eigen::Dynamic, 2> total_cs_matrix_;
164 Interpolator::LogLogLinear total_cs_interpolator_;
165 Eigen::Matrix<double, Eigen::Dynamic, 2> incoherent_scattering_function_matrix_;
166 Interpolator::LogLogLinear incoherent_scattering_function_interpolator_;
167 Eigen::Matrix<double, Eigen::Dynamic, 2> coherent_form_factor_matrix_;
168 Interpolator::LogLogLinear coherent_form_factor_interpolator_;
169 Eigen::Matrix<double, Eigen::Dynamic, 2> mass_energy_absorption_coefficient_matrix_;
170 Interpolator::LogLogLinear mass_energy_absorption_coefficient_interpolator_;
171 ProbabilityDist::ContinuousInversion coherent_scattering_dcs_dist_;
172
173 void initializeData();
174
175 void setInteractionCrossSectionsAndInterpolators();
176 void setTotalCrossSectionsAndInterpolator();
177 void setIncoherentScatteringFunctionAndInterpolator();
178 void setCoherentScatteringFormFactorAndInterpolator();
179 void setMassEnergyAbsorptionCoefficientsAndInterpolator();
180
181 void setIncoherentScatteringCrossSectionAndInterpolator();
182 void setCoherentScatteringCrossSectionAndInterpolator();
183 void setPhotoelectricCrossSectionAndInterpolator();
184
185 void setCoherentScatteringDCSDistribution();
186
187 Eigen::Matrix<double, Eigen::Dynamic, 2> getTotalCrossSectionsMatrixFromInteractionData();
188 Eigen::Matrix<double, Eigen::Dynamic, 2> calculateWeightedAverageOfColumns(const std::string &tableName, const std::string &dataColumnName,
189 bool scale_to_macroscopic = false);
190
191 void fillTotalCrossSectionsMatrix(Eigen::MatrixXd& total_cross_sections_matrix, const Eigen::MatrixXd& merged_energy_matrix);
192 std::unordered_map<int, Eigen::Matrix<double, Eigen::Dynamic, 2>> getTableMatrixForAllElements(const std::string &tableName, const std::string &dataColumnName);
193 std::unordered_map<int, std::shared_ptr<Interpolator::Interpolator>> getInterpolatorsForAllElements(const std::string &tableName,
194 const std::unordered_map<int, Eigen::Matrix<double, Eigen::Dynamic, 2>>& table_matrix_map);
195
196 static std::shared_ptr<Interpolator::Interpolator> getInterpolatorForElement(const std::string& tableName, const Eigen::Matrix<double, Eigen::Dynamic, 2>& matrix);
197
198 Eigen::Matrix<double, Eigen::Dynamic, 2> getTableMatrix(const std::string& tableName, const std::string& dataColumnName, int element);
199};
200
201#endif //MCXRAYTRANSPORT_MATERIAL_DATA_H
Class which provides an interface to a SQLite database.
Class which performs linear interpolation on a log-log scale.
Class which performs 3rd order spline interpolation on a log-log scale.
Class which contains the computed simulation data for a material.
Eigen::Matrix< double, Eigen::Dynamic, 2 > getCoherentScatteringCrossSectionMatrix() const
Returns the 2 column matrix containing the coherent scattering cross sections for the material.
double sampleCoherentScatteringDCS(double energy) const
Samples the incoherent scattering DCS distribution for the material at the given energy.
double interpolateMassEnergyAbsorptionCoefficient(double energy) const
Returns the interpolated mass energy absorption coefficient for the material at the given energy.
Eigen::Matrix< double, Eigen::Dynamic, 2 > getPhotoelectricCrossSectionMatrix() const
Returns the 2 column matrix containing the photoelectric cross sections for the material.
double interpolateCoherentFormFactor(double x) const
Returns the interpolated coherent scattering form factor for the material at the given momentum trans...
MaterialData(MaterialProperties &properties, DataAccessObject &dao)
Constructor for the MaterialData class.
Eigen::Matrix< double, Eigen::Dynamic, 2 > getTotalCrossSectionMatrix() const
Returns the 2 column matrix containing the total cross sections for the material.
double interpolateTotalCrossSection(double energy) const
Returns the interpolated total cross section for the material at the given energy.
double interpolateCoherentScatteringCrossSection(double energy) const
Returns the interpolated coherent scattering cross section for the material at the given energy.
double interpolateIncoherentScatteringFunction(double x) const
Returns the interpolated incoherent scattering function for the material at the given momentum transf...
double interpolatePhotoelectricCrossSection(double energy) const
Returns the interpolated photoelectric cross section for the material at the given energy.
Eigen::Matrix< double, Eigen::Dynamic, 2 > getCoherentFormFactorMatrix() const
Returns the 2 column matrix containing the coherent scattering form factors for the material.
Eigen::Matrix< double, Eigen::Dynamic, 2 > getIncoherentScatteringCrossSectionMatrix() const
Returns the 2 column matrix containing the incoherent scattering cross sections for the material.
Eigen::Matrix< double, Eigen::Dynamic, 2 > getIncoherentScatteringFunctionMatrix() const
Returns the 2 column matrix containing the incoherent scattering functions for the material.
double interpolateIncoherentScatteringCrossSection(double energy) const
Returns the interpolated incoherent scattering cross section for the material at the given energy.
Class which represents the properties of a material.
Class which uses PENELOPE's RITA inversion sampling algorithm to sample from a continuous distributio...
double sample(double E) const
Returns a sample from the continuous distribution.
Contains constants used throughout the code.