MIDSX 0.1
A x-ray transport code system for dosimetry
Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
MaterialData Class Reference

Class which contains the computed simulation data for a material. More...

#include <material_data.h>

Collaboration diagram for MaterialData:

Public Member Functions

 MaterialData (MaterialProperties &properties, DataAccessObject &dao)
 Constructor for the MaterialData class.
 
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 > getCoherentScatteringCrossSectionMatrix () const
 Returns the 2 column matrix containing the coherent scattering cross sections for the material.
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > getPhotoelectricCrossSectionMatrix () const
 Returns the 2 column matrix containing the photoelectric cross sections for the material.
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > getTotalCrossSectionMatrix () const
 Returns the 2 column matrix containing the total 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.
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > getCoherentFormFactorMatrix () const
 Returns the 2 column matrix containing the coherent scattering form factors for the material.
 
double interpolateIncoherentScatteringCrossSection (double energy) const
 Returns the interpolated incoherent scattering 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 interpolatePhotoelectricCrossSection (double energy) const
 Returns the interpolated photoelectric cross section for the material at the given energy.
 
double interpolateTotalCrossSection (double energy) const
 Returns the interpolated total 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 transfer.
 
double interpolateCoherentFormFactor (double x) const
 Returns the interpolated coherent scattering form factor for the material at the given momentum transfer.
 
double interpolateMassEnergyAbsorptionCoefficient (double energy) const
 Returns the interpolated mass energy absorption coefficient for the material at the given energy.
 
double sampleCoherentScatteringDCS (double energy) const
 Samples the incoherent scattering DCS distribution for the material at the given energy.
 

Private Member Functions

void initializeData ()
 
void setInteractionCrossSectionsAndInterpolators ()
 
void setTotalCrossSectionsAndInterpolator ()
 
void setIncoherentScatteringFunctionAndInterpolator ()
 
void setCoherentScatteringFormFactorAndInterpolator ()
 
void setMassEnergyAbsorptionCoefficientsAndInterpolator ()
 
void setIncoherentScatteringCrossSectionAndInterpolator ()
 
void setCoherentScatteringCrossSectionAndInterpolator ()
 
void setPhotoelectricCrossSectionAndInterpolator ()
 
void setCoherentScatteringDCSDistribution ()
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > getTotalCrossSectionsMatrixFromInteractionData ()
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > calculateWeightedAverageOfColumns (const std::string &tableName, const std::string &dataColumnName, bool scale_to_macroscopic=false)
 
void fillTotalCrossSectionsMatrix (Eigen::MatrixXd &total_cross_sections_matrix, const Eigen::MatrixXd &merged_energy_matrix)
 
std::unordered_map< int, Eigen::Matrix< double, Eigen::Dynamic, 2 > > getTableMatrixForAllElements (const std::string &tableName, const std::string &dataColumnName)
 
std::unordered_map< int, std::shared_ptr< Interpolator::Interpolator > > getInterpolatorsForAllElements (const std::string &tableName, const std::unordered_map< int, Eigen::Matrix< double, Eigen::Dynamic, 2 > > &table_matrix_map)
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > getTableMatrix (const std::string &tableName, const std::string &dataColumnName, int element)
 

Static Private Member Functions

static std::shared_ptr< Interpolator::InterpolatorgetInterpolatorForElement (const std::string &tableName, const Eigen::Matrix< double, Eigen::Dynamic, 2 > &matrix)
 

Private Attributes

MaterialPropertiesproperties_
 
DataAccessObjectdao_
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > incoherent_cs_matrix_
 
Interpolator::LogLogSpline incoherent_cs_interpolator_
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > coherent_cs_matrix_
 
Interpolator::LogLogSpline coherent_cs_interpolator_
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > photoelectric_cs_matrix_
 
Interpolator::LogLogLinear photoelectric_cs_interpolator_
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > total_cs_matrix_
 
Interpolator::LogLogLinear total_cs_interpolator_
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > incoherent_scattering_function_matrix_
 
Interpolator::LogLogLinear incoherent_scattering_function_interpolator_
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > coherent_form_factor_matrix_
 
Interpolator::LogLogLinear coherent_form_factor_interpolator_
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > mass_energy_absorption_coefficient_matrix_
 
Interpolator::LogLogLinear mass_energy_absorption_coefficient_interpolator_
 
ProbabilityDist::ContinuousInversion coherent_scattering_dcs_dist_
 

Detailed Description

Class which contains the computed simulation data for a material.

This class contains the computed data for a material obtained from the SQLite database.

Definition at line 19 of file material_data.h.

Constructor & Destructor Documentation

◆ MaterialData()

MaterialData::MaterialData ( MaterialProperties & properties,
DataAccessObject & dao )

Constructor for the MaterialData class.

Parameters
propertiesThe MaterialProperties object for the material.
daoThe DataAccessObject to use to query the database.

Member Function Documentation

◆ getCoherentFormFactorMatrix()

Eigen::Matrix< double, Eigen::Dynamic, 2 > MaterialData::getCoherentFormFactorMatrix ( ) const
inline

Returns the 2 column matrix containing the coherent scattering form factors for the material.

The first column contains the momentum transfer values, and the second column contains the coherent scattering form factor. Momentum transfer values are in inverse angstroms and the form factor is unitless.

Returns
The 2 column matrix containing the coherent scattering form factors for the material.

Definition at line 87 of file material_data.h.

◆ getCoherentScatteringCrossSectionMatrix()

Eigen::Matrix< double, Eigen::Dynamic, 2 > MaterialData::getCoherentScatteringCrossSectionMatrix ( ) const
inline

Returns the 2 column matrix containing the coherent scattering cross sections for the material.

The first column contains the energy values, and the second column contains the coherent scattering cross sections. Energy values are in eV and cross sections are in cm^-1.

Returns
The 2 column matrix containing the coherent scattering cross sections for the material.

Definition at line 47 of file material_data.h.

◆ getIncoherentScatteringCrossSectionMatrix()

Eigen::Matrix< double, Eigen::Dynamic, 2 > MaterialData::getIncoherentScatteringCrossSectionMatrix ( ) const
inline

Returns the 2 column matrix containing the incoherent scattering cross sections for the material.

The first column contains the energy values, and the second column contains the incoherent scattering cross sections. Energy values are in eV and cross sections are in cm^-1.

Returns
The 2 column matrix containing the incoherent scattering cross sections for the material.

Definition at line 37 of file material_data.h.

◆ getIncoherentScatteringFunctionMatrix()

Eigen::Matrix< double, Eigen::Dynamic, 2 > MaterialData::getIncoherentScatteringFunctionMatrix ( ) const
inline

Returns the 2 column matrix containing the incoherent scattering functions for the material.

The first column contains the momentum transfer values, and the second column contains the incoherent scattering function. Momentum transfer values are in inverse angstroms and the scattering function is unitless.

Returns
The 2 column matrix containing the incoherent scattering functions for the material.

Definition at line 77 of file material_data.h.

◆ getPhotoelectricCrossSectionMatrix()

Eigen::Matrix< double, Eigen::Dynamic, 2 > MaterialData::getPhotoelectricCrossSectionMatrix ( ) const
inline

Returns the 2 column matrix containing the photoelectric cross sections for the material.

The first column contains the energy values, and the second column contains the photoelectric cross sections. Energy values are in eV and cross sections are in cm^-1.

Returns
The 2 column matrix containing the photoelectric cross sections for the material.

Definition at line 57 of file material_data.h.

◆ getTotalCrossSectionMatrix()

Eigen::Matrix< double, Eigen::Dynamic, 2 > MaterialData::getTotalCrossSectionMatrix ( ) const
inline

Returns the 2 column matrix containing the total cross sections for the material.

The first column contains the energy values, and the second column contains the total cross sections. Energy values are in eV and cross sections are in cm^-1.

Returns
The 2 column matrix containing the total cross sections for the material.

Definition at line 67 of file material_data.h.

◆ interpolateCoherentFormFactor()

double MaterialData::interpolateCoherentFormFactor ( double x) const
inline

Returns the interpolated coherent scattering form factor for the material at the given momentum transfer.

Parameters
xThe momentum transfer (inverse angstroms) at which to return the coherent scattering form factor.
Returns
The interpolated coherent scattering form factor (unitless) for the material at the given momentum transfer.

Definition at line 135 of file material_data.h.

◆ interpolateCoherentScatteringCrossSection()

double MaterialData::interpolateCoherentScatteringCrossSection ( double energy) const
inline

Returns the interpolated coherent scattering cross section for the material at the given energy.

Parameters
energyThe energy (eV) at which to return the coherent scattering cross section.
Returns
The interpolated coherent scattering cross section (cm^-1) for the material at the given energy.

Definition at line 103 of file material_data.h.

◆ interpolateIncoherentScatteringCrossSection()

double MaterialData::interpolateIncoherentScatteringCrossSection ( double energy) const
inline

Returns the interpolated incoherent scattering cross section for the material at the given energy.

Parameters
energyThe energy (eV) at which to return the incoherent scattering cross section.
Returns
The interpolated incoherent scattering cross section (cm^-1) for the material at the given energy.

Definition at line 95 of file material_data.h.

◆ interpolateIncoherentScatteringFunction()

double MaterialData::interpolateIncoherentScatteringFunction ( double x) const
inline

Returns the interpolated incoherent scattering function for the material at the given momentum transfer.

Parameters
xThe momentum transfer (inverse angstroms) at which to return the incoherent scattering function.
Returns
The interpolated incoherent scattering function (unitless) for the material at the given momentum transfer.

Definition at line 127 of file material_data.h.

◆ interpolateMassEnergyAbsorptionCoefficient()

double MaterialData::interpolateMassEnergyAbsorptionCoefficient ( double energy) const
inline

Returns the interpolated mass energy absorption coefficient for the material at the given energy.

Parameters
energyThe energy (eV) at which to return the mass energy absorption coefficient.
Returns
The interpolated mass energy absorption coefficient (cm^2/g) for the material at the given energy.

Definition at line 143 of file material_data.h.

◆ interpolatePhotoelectricCrossSection()

double MaterialData::interpolatePhotoelectricCrossSection ( double energy) const
inline

Returns the interpolated photoelectric cross section for the material at the given energy.

Parameters
energyThe energy (eV) at which to return the photoelectric cross section.
Returns
The interpolated photoelectric cross section (cm^-1) for the material at the given energy.

Definition at line 111 of file material_data.h.

◆ interpolateTotalCrossSection()

double MaterialData::interpolateTotalCrossSection ( double energy) const
inline

Returns the interpolated total cross section for the material at the given energy.

Parameters
energyThe energy (eV) at which to return the total cross section.
Returns
The interpolated total cross section (cm^-1) for the material at the given energy.

Definition at line 119 of file material_data.h.

◆ sampleCoherentScatteringDCS()

double MaterialData::sampleCoherentScatteringDCS ( double energy) const
inline

Samples the incoherent scattering DCS distribution for the material at the given energy.

Parameters
energyThe energy (eV) at which to sample the incoherent scattering DCS distribution.
Returns
The sampled incoherent scattering DCS for the material at the given energy.

Definition at line 151 of file material_data.h.

Here is the call graph for this function:

Member Data Documentation

◆ coherent_cs_interpolator_

Interpolator::LogLogSpline MaterialData::coherent_cs_interpolator_
private

Definition at line 160 of file material_data.h.

◆ coherent_cs_matrix_

Eigen::Matrix<double, Eigen::Dynamic, 2> MaterialData::coherent_cs_matrix_
private

Definition at line 159 of file material_data.h.

◆ coherent_form_factor_interpolator_

Interpolator::LogLogLinear MaterialData::coherent_form_factor_interpolator_
private

Definition at line 168 of file material_data.h.

◆ coherent_form_factor_matrix_

Eigen::Matrix<double, Eigen::Dynamic, 2> MaterialData::coherent_form_factor_matrix_
private

Definition at line 167 of file material_data.h.

◆ coherent_scattering_dcs_dist_

ProbabilityDist::ContinuousInversion MaterialData::coherent_scattering_dcs_dist_
private

Definition at line 171 of file material_data.h.

◆ dao_

DataAccessObject& MaterialData::dao_
private

Definition at line 155 of file material_data.h.

◆ incoherent_cs_interpolator_

Interpolator::LogLogSpline MaterialData::incoherent_cs_interpolator_
private

Definition at line 158 of file material_data.h.

◆ incoherent_cs_matrix_

Eigen::Matrix<double, Eigen::Dynamic, 2> MaterialData::incoherent_cs_matrix_
private

Definition at line 157 of file material_data.h.

◆ incoherent_scattering_function_interpolator_

Interpolator::LogLogLinear MaterialData::incoherent_scattering_function_interpolator_
private

Definition at line 166 of file material_data.h.

◆ incoherent_scattering_function_matrix_

Eigen::Matrix<double, Eigen::Dynamic, 2> MaterialData::incoherent_scattering_function_matrix_
private

Definition at line 165 of file material_data.h.

◆ mass_energy_absorption_coefficient_interpolator_

Interpolator::LogLogLinear MaterialData::mass_energy_absorption_coefficient_interpolator_
private

Definition at line 170 of file material_data.h.

◆ mass_energy_absorption_coefficient_matrix_

Eigen::Matrix<double, Eigen::Dynamic, 2> MaterialData::mass_energy_absorption_coefficient_matrix_
private

Definition at line 169 of file material_data.h.

◆ photoelectric_cs_interpolator_

Interpolator::LogLogLinear MaterialData::photoelectric_cs_interpolator_
private

Definition at line 162 of file material_data.h.

◆ photoelectric_cs_matrix_

Eigen::Matrix<double, Eigen::Dynamic, 2> MaterialData::photoelectric_cs_matrix_
private

Definition at line 161 of file material_data.h.

◆ properties_

MaterialProperties& MaterialData::properties_
private

Definition at line 154 of file material_data.h.

◆ total_cs_interpolator_

Interpolator::LogLogLinear MaterialData::total_cs_interpolator_
private

Definition at line 164 of file material_data.h.

◆ total_cs_matrix_

Eigen::Matrix<double, Eigen::Dynamic, 2> MaterialData::total_cs_matrix_
private

Definition at line 163 of file material_data.h.


The documentation for this class was generated from the following file: