CTL  0.6.1
Computed Tomography Library
Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
CTL::SpectralVolumeData Class Reference

The SpectralVolumeData class holds voxelized data (for a single material) along with information on its spectrally-dependent absorption properties. More...

#include <spectralvolumedata.h>

Inheritance diagram for CTL::SpectralVolumeData:
Inheritance graph
[legend]
Collaboration diagram for CTL::SpectralVolumeData:
Collaboration graph
[legend]

Public Member Functions

 SpectralVolumeData (VoxelVolume< float > muValues)
 Constructs a SpectralVolumeData representing the attenuation coefficients muValues (in 1/mm). More...
 
 SpectralVolumeData (VoxelVolume< float > muValues, std::shared_ptr< AbstractIntegrableDataModel > absorptionModel, float referenceEnergy, const QString &materialName=QString())
 Constructs a SpectralVolumeData representing the attenuation coefficients muValues (in 1/mm) with respect to the given referenceEnergy (in keV) for a material described by absorptionModel. More...
 
 SpectralVolumeData (VoxelVolume< float > materialDensity, std::shared_ptr< AbstractIntegrableDataModel > absorptionModel, const QString &materialName=QString())
 Constructs a SpectralVolumeData representing the density values materialDensity (in g/cm^3) for a material described by absorptionModel. More...
 
virtual SpectralVolumeDataclone () const
 
std::shared_ptr< AbstractIntegrableDataModelabsorptionModel () const
 
std::unique_ptr< SpectralVolumeDatadensityVolume () const
 
bool hasSpectralInformation () const
 
bool isDensityVolume () const
 
bool isMuVolume () const
 
float massAttenuationCoeff (float atEnergy) const
 
const QString & materialName () const
 
float meanMassAttenuationCoeff (float centerEnergy, float binWidth) const
 
std::unique_ptr< SpectralVolumeDatamuVolume (float referenceEnergy) const
 
std::unique_ptr< SpectralVolumeDatamuVolume (float centerEnergy, float binWidth) const
 
float referenceEnergy () const
 
float referenceMassAttenuationCoeff () const
 
void replaceAbsorptionModel (AbstractIntegrableDataModel *absorptionModel)
 
void replaceAbsorptionModel (std::shared_ptr< AbstractIntegrableDataModel > absorptionModel)
 
void setDensity (VoxelVolume< float > density)
 
void setMaterialName (const QString &name)
 
 SpectralVolumeData (const SpectralVolumeData &)=default
 
 SpectralVolumeData (SpectralVolumeData &&)=default
 
SpectralVolumeDataoperator= (const SpectralVolumeData &)=default
 
SpectralVolumeDataoperator= (SpectralVolumeData &&)=default
 
void transformToAttenuationCoeff (float referenceEnergy)
 
void transformToDensity ()
 
- Public Member Functions inherited from CTL::VoxelVolume< float >
iterator begin ()
 
const_iterator begin () const
 
iterator end ()
 
const_iterator end () const
 
const_iterator cbegin () const
 
const_iterator cend () const
 
reverse_iterator rbegin ()
 
const_reverse_iterator rbegin () const
 
reverse_iterator rend ()
 
const_reverse_iterator rend () const
 
const_reverse_iterator crbegin () const
 
const_reverse_iterator crend () const
 
 VoxelVolume (const Dimensions &nbVoxels)
 
 VoxelVolume (const Dimensions &nbVoxels, const VoxelSize &voxelSize)
 
 VoxelVolume (uint nbVoxelX, uint nbVoxelY, uint nbVoxelZ)
 
 VoxelVolume (uint nbVoxelX, uint nbVoxelY, uint nbVoxelZ, float xSize, float ySize, float zSize)
 
 VoxelVolume (const Dimensions &nbVoxels, std::vector< float > data)
 
 VoxelVolume (const Dimensions &nbVoxels, const VoxelSize &voxelSize, std::vector< float > data)
 
 VoxelVolume (uint nbVoxelX, uint nbVoxelY, uint nbVoxelZ, std::vector< float > data)
 
 VoxelVolume (uint nbVoxelX, uint nbVoxelY, uint nbVoxelZ, float xSize, float ySize, float zSize, std::vector< float > data)
 
 VoxelVolume (const VoxelVolume &)=default
 
 VoxelVolume (VoxelVolume &&)=default
 
VoxelVolumeoperator= (const VoxelVolume &)=default
 
VoxelVolumeoperator= (VoxelVolume &&)=default
 
size_t allocatedElements () const
 
const std::vector< float > & constData () const
 
const std::vector< float > & data () const
 
std::vector< float > & data ()
 
const Dimensionsdimensions () const
 
bool hasData () const
 
const DimensionsnbVoxels () const
 
const Offsetoffset () const
 
float * rawData ()
 
const float * rawData () const
 
size_t totalVoxelCount () const
 
const VoxelSizevoxelSize () const
 
void setData (std::vector< float > &&data)
 
void setData (const std::vector< float > &data)
 
void setVolumeOffset (const Offset &offset)
 
void setVolumeOffset (float xOffset, float yOffset, float zOffset)
 
void setVoxelSize (const VoxelSize &size)
 
void setVoxelSize (float xSize, float ySize, float zSize)
 
void setVoxelSize (float isotropicSize)
 
void allocateMemory ()
 
void allocateMemory (const float &initValue)
 
VoxelCoordinates coordinates (const VoxelIndex &index) const
 
VoxelCoordinates coordinates (uint x, uint y, uint z) const
 
VoxelCoordinates cornerVoxel () const
 
VoxelIndex index (const VoxelCoordinates &coordinates) const
 
VoxelIndex index (float x_mm, float y_mm, float z_mm) const
 
bool isIsotropic () const
 
void fill (const float &fillValue)
 
void freeMemory ()
 
float max () const
 
float min () const
 
VoxelVolume< float > reslicedByX (bool reverse=false) const
 
VoxelVolume< float > reslicedByY (bool reverse=false) const
 
VoxelVolume< float > reslicedByZ (bool reverse=false) const
 
Chunk2D< float > sliceX (uint slice) const
 
Chunk2D< float > sliceY (uint slice) const
 
Chunk2D< float > sliceZ (uint slice) const
 
float smallestVoxelSize () const
 
VoxelCoordinates volumeCorner () const
 
std::vector< float >::reference operator() (uint x, uint y, uint z)
 
std::vector< float >::const_reference operator() (uint x, uint y, uint z) const
 
std::vector< float >::reference operator() (const VoxelIndex &index)
 
std::vector< float >::const_reference operator() (const VoxelIndex &index) const
 
VoxelVolume< float > & operator+= (const VoxelVolume< float > &other)
 
VoxelVolume< float > & operator+= (const float &additiveShift)
 
VoxelVolume< float > & operator-= (const VoxelVolume< float > &other)
 
VoxelVolume< float > & operator-= (const float &subtractiveShift)
 
VoxelVolume< float > & operator *= (const float &factor)
 
VoxelVolume< float > & operator/= (const float &divisor)
 
VoxelVolume< float > operator+ (const VoxelVolume< float > &other) const
 
VoxelVolume< float > operator+ (const float &additiveShift) const
 
VoxelVolume< float > operator- (const VoxelVolume< float > &other) const
 
VoxelVolume< float > operator- (const float &subtractiveShift) const
 
VoxelVolume< float > operator * (const float &factor) const
 
VoxelVolume< float > operator/ (const float &divisor) const
 

Static Public Member Functions

static SpectralVolumeData ball (float radius, float voxelSize, float density, std::shared_ptr< AbstractIntegrableDataModel > absorptionModel)
 
static SpectralVolumeData cube (uint nbVoxel, float voxelSize, float density, std::shared_ptr< AbstractIntegrableDataModel > absorptionModel)
 
static SpectralVolumeData cylinderX (float radius, float height, float voxelSize, float density, std::shared_ptr< AbstractIntegrableDataModel > absorptionModel)
 
static SpectralVolumeData cylinderY (float radius, float height, float voxelSize, float density, std::shared_ptr< AbstractIntegrableDataModel > absorptionModel)
 
static SpectralVolumeData cylinderZ (float radius, float height, float voxelSize, float density, std::shared_ptr< AbstractIntegrableDataModel > absorptionModel)
 
static SpectralVolumeData fromMuVolume (VoxelVolume< float > muValues, std::shared_ptr< AbstractIntegrableDataModel > absorptionModel, float referenceEnergy=50.0f)
 
static SpectralVolumeData fromHUVolume (VoxelVolume< float > HUValues, std::shared_ptr< AbstractIntegrableDataModel > absorptionModel, float referenceEnergy=50.0f)
 
- Static Public Member Functions inherited from CTL::VoxelVolume< float >
static VoxelVolume< float > fromChunk2DStack (const std::vector< Chunk2D< float >> &stack)
 
static VoxelVolume< float > ball (float radius, float voxelSize, const float &fillValue)
 
static VoxelVolume< float > cube (uint nbVoxel, float voxelSize, const float &fillValue)
 
static VoxelVolume< float > cylinderX (float radius, float height, float voxelSize, const float &fillValue)
 
static VoxelVolume< float > cylinderY (float radius, float height, float voxelSize, const float &fillValue)
 
static VoxelVolume< float > cylinderZ (float radius, float height, float voxelSize, const float &fillValue)
 

Private Member Functions

void changeReferenceEnergy (float newReferenceEnergy)
 
void changeReferenceMassAttCoeff (float newReferenceMassAttCoeff, float correspondingRefEnergy)
 
void transformToAttenuationCoeff (float referenceMassAttCoeff, float correspondingRefEnergy)
 

Private Attributes

std::shared_ptr< AbstractIntegrableDataModel_absorptionModel
 
QString _materialName
 
bool _hasNonDefaultAbsModel { false }
 True if absorption model has been set explicitely.
 
bool _isMu { false }
 True if voxel data represents att. coefficents.
 
float _refEnergy { -1.0f }
 Reference energy corresponding to mu values.
 
float _refMassAttenuationCoeff { -1.0f }
 Reference att. coeff. corresponding to mu values.
 

Additional Inherited Members

- Public Types inherited from CTL::VoxelVolume< float >
using Dimensions = VoxelVolumeDimensions
 
using VoxelSize = VoxelVolumeVoxelSize
 
using Offset = VoxelVolumeOffset
 
using iterator = VoxelIterator< typename std::vector< float >::iterator >
 
using const_iterator = VoxelIterator< typename std::vector< float >::const_iterator >
 
using reverse_iterator = std::reverse_iterator< iterator >
 
using const_reverse_iterator = std::reverse_iterator< const_iterator >
 
- Protected Attributes inherited from CTL::VoxelVolume< float >
Dimensions _dim
 The dimensions of the volume.
 
VoxelSize _size
 The size of individual voxels (in mm).
 
Offset _offset
 The positional offset of the volume (in mm).
 
std::vector< float > _data
 The internal data of the volume.
 

Detailed Description

The SpectralVolumeData class holds voxelized data (for a single material) along with information on its spectrally-dependent absorption properties.

This class is an extension of the VoxelVolume class (template type float) by means of enhancing it with information on its spectrally-dependent absorption properties. It represents a single material (either elemental or composite), in the way that the spectral dependency of the mass attenuation coefficient is assumed identical for all voxels (use CompositeVolume in cases where this assumption is too restrictive).

Spectral information always requires a data model describing the spectral dependency of the mass-attenuation coefficient for the material represented by this volume. Additionally, one of the following conditions must be fulfilled:

The two different representations can be queried explicitely with densityVolume() and muVolume(). Note that both these methods need to create copies of the data.

Availability of spectral information in the object can be checked via hasSpectralInformation(). In case the abovementioned conditions are not entirely fulfilled, the stored data is considered to represent straightforward attenuation coefficients (mu, in 1/mm) without any spectral information. Internally, such a volume object works with a constant attenuation model (i.e. mu/rho = 1.0 cm^2/g across entire energy range). This behavior is implemented only for convenience and it is recommended to use plain VoxelVolume<float> objects for such cases. Any VoxelVolume<float> object representing attenuation values (either in 1/mm (i.e. mu) or Hounsfield units) can be transformed into SpectralVolumeData using the factory methods fromMuVolume() or fromHUVolume() to complement the missing spectral information.

Constructor & Destructor Documentation

◆ SpectralVolumeData() [1/3]

CTL::SpectralVolumeData::SpectralVolumeData ( VoxelVolume< float >  muValues)

Constructs a SpectralVolumeData representing the attenuation coefficients muValues (in 1/mm).

This is a convenience constructor, mainly intended to allow for implicit casts of VoxelVolume<float> to SpectralVolumeData.

Note that it is not meaningful to use this constructor with input data in Hounsfield units (HU).

No spectral information is available in the resulting object. It is strongly encouraged to use VoxelVolume<float> directly when managing non-spectral attenuation information.

◆ SpectralVolumeData() [2/3]

CTL::SpectralVolumeData::SpectralVolumeData ( VoxelVolume< float >  muValues,
std::shared_ptr< AbstractIntegrableDataModel absorptionModel,
float  referenceEnergy,
const QString &  materialName = QString() 
)

Constructs a SpectralVolumeData representing the attenuation coefficients muValues (in 1/mm) with respect to the given referenceEnergy (in keV) for a material described by absorptionModel.

Note that this can only be used when muValues holds values in 1/mm. To create from Hounsfield units (HU), use fromHUVolume().

Keeps internally stored data in attenuation domain. The factory method fromMuVolume() has similar function, but transforms data to density domain.

◆ SpectralVolumeData() [3/3]

CTL::SpectralVolumeData::SpectralVolumeData ( VoxelVolume< float >  materialDensity,
std::shared_ptr< AbstractIntegrableDataModel absorptionModel,
const QString &  materialName = QString() 
)

Constructs a SpectralVolumeData representing the density values materialDensity (in g/cm^3) for a material described by absorptionModel.

Member Function Documentation

◆ absorptionModel()

std::shared_ptr< AbstractIntegrableDataModel > CTL::SpectralVolumeData::absorptionModel ( ) const

Returns the absorption model of this instance.

The absorption model represents a single material (either elemental or composite) and contains the spectral dependency of its mass attenuation coefficients (in cm^2/g).

◆ ball()

SpectralVolumeData CTL::SpectralVolumeData::ball ( float  radius,
float  voxelSize,
float  density,
std::shared_ptr< AbstractIntegrableDataModel absorptionModel 
)
static

Creates a SpectralVolumeData object that represents a voxelized ball with radius radius, isometric voxel size voxelSize (both in mm) and filled (homogeneuosly) with density value density (in g/cm^3). The material properties (i.e. spectrally-dependent mass attenuation coefficients) are specified by absorptionModel. The voxels surrounding the ball are filled with density 0.0 g/cm^3.

◆ cube()

SpectralVolumeData CTL::SpectralVolumeData::cube ( uint  nbVoxel,
float  voxelSize,
float  density,
std::shared_ptr< AbstractIntegrableDataModel absorptionModel 
)
static

Constructs a cubic SpectralVolumeData with nbVoxel x nbVoxel x nbVoxel voxels (voxel dimension: voxelSize x voxelSize x voxelSize), filled (homogeneuosly) with density value density (in g/cm^3). The material properties (i.e. spectrally-dependent mass attenuation coefficients) are specified by absorptionModel.

◆ cylinderX()

SpectralVolumeData CTL::SpectralVolumeData::cylinderX ( float  radius,
float  height,
float  voxelSize,
float  density,
std::shared_ptr< AbstractIntegrableDataModel absorptionModel 
)
static

Creates a SpectralVolumeData object that represents a voxelized cylinder with radius radius and height height (both in mm) that is aligned with the x-axis. It has isometric voxel size voxelSize (in mm) and is filled (homogeneuosly) with density value density (in g/cm^3). The material properties (i.e. spectrally-dependent mass attenuation coefficients) are specified by absorptionModel. The voxels surrounding the cylinder are filled with density 0.0 g/cm^3.

The resulting volume will have \( \left\lceil 2\cdot radius/voxelSize\right\rceil \) voxels in y- and z-dimension and \( \left\lceil height/voxelSize\right\rceil \) in x-direction.

◆ cylinderY()

SpectralVolumeData CTL::SpectralVolumeData::cylinderY ( float  radius,
float  height,
float  voxelSize,
float  density,
std::shared_ptr< AbstractIntegrableDataModel absorptionModel 
)
static

Creates a SpectralVolumeData object that represents a voxelized cylinder with radius radius and height height (both in mm) that is aligned with the y-axis. It has isometric voxel size voxelSize (in mm) and is filled (homogeneuosly) with density value density (in g/cm^3). The material properties (i.e. spectrally-dependent mass attenuation coefficients) are specified by absorptionModel. The voxels surrounding the cylinder are filled with density 0.0 g/cm^3.

The resulting volume will have \( \left\lceil 2\cdot radius/voxelSize\right\rceil \) voxels in x- and z-dimension and \( \left\lceil height/voxelSize\right\rceil \) in y-direction.

◆ cylinderZ()

SpectralVolumeData CTL::SpectralVolumeData::cylinderZ ( float  radius,
float  height,
float  voxelSize,
float  density,
std::shared_ptr< AbstractIntegrableDataModel absorptionModel 
)
static

Creates a SpectralVolumeData object that represents a voxelized cylinder with radius radius and height height (both in mm) that is aligned with the z-axis. It has isometric voxel size voxelSize (in mm) and is filled (homogeneuosly) with density value density (in g/cm^3). The material properties (i.e. spectrally-dependent mass attenuation coefficients) are specified by absorptionModel. The voxels surrounding the cylinder are filled with density 0.0 g/cm^3.

The resulting volume will have \( \left\lceil 2\cdot radius/voxelSize\right\rceil \) voxels in x- and y-dimension and \( \left\lceil height/voxelSize\right\rceil \) in z-direction.

◆ densityVolume()

std::unique_ptr< SpectralVolumeData > CTL::SpectralVolumeData::densityVolume ( ) const

Returns the density representation of this instance.

Note that this creates a copy of the data. In case this instance does already contain density data (check this via isDensityVolume()) it is discouraged to call this method because it would only create an unnecessary copy.

◆ fromHUVolume()

SpectralVolumeData CTL::SpectralVolumeData::fromHUVolume ( VoxelVolume< float >  HUValues,
std::shared_ptr< AbstractIntegrableDataModel absorptionModel,
float  referenceEnergy = 50.0f 
)
static

Creates a SpectralVolumeData object from the attenuation values given by HUValues (in Hounsfield units) representing the material specified by its spectrally-dependent mass attenuation coefficients in absorptionModel. For meaningful results, you need to also specify the reference energy, to which the Hounsfield units correspond, by referenceEnergy (in keV).

Generates the density representation of the data.

// Starting point: voxel data in Hounsfield units (eg. loaded patient data)
// here: make a cube phantom (100^3 voxels, size 1 mm^3) filled with value 100 HU
auto volumeHU = VoxelVolume<float>(100, 100, 100, 1.0f, 1.0f, 1.0f);
volumeHU.fill(100.0f); // fill with value 100 HU
// Now: create spectral volume from HU data (soft tissue data, reference energy of 62 keV)
auto volume = SpectralVolumeData::fromHUVolume(volumeHU, /* std::move(volumeHU), if no longer required. */
database::attenuationModel(database::Composite::Tissue_Soft),
62.0f);
qInfo() << volume.max(); // output: 1.10614 (g/cm^3) <-- density representation

◆ fromMuVolume()

SpectralVolumeData CTL::SpectralVolumeData::fromMuVolume ( VoxelVolume< float >  muValues,
std::shared_ptr< AbstractIntegrableDataModel absorptionModel,
float  referenceEnergy = 50.0f 
)
static

Creates a SpectralVolumeData object from the attenuation values given by muValues (in 1 / mm) corresponding to the reference energy referenceEnergy of the material specified by its spectrally-dependent mass attenuation coefficients in absorptionModel.

Generates the density representation of the data. To prevent transformation into density domain (e.g. if follow-up processing needs to be done in attenuation domain anyway), use the constructor with the same input instead.

// Starting point: voxel data with attenuation coefficients in 1/mm
// here: make a cube phantom (100^3 voxels, size 1 mm^3) representing cortical bone
auto volumeMu = VoxelVolume<float>(100, 100, 100, 1.0f, 1.0f, 1.0f);
volumeMu.fill(0.056f); // fill with value 0.056/mm <-- approx. attenuation coeff. of bone @ 65 keV
// create spectral volume from mu (i.e. att. coeff.) data (cortical bone data, reference energy of 65 keV)
auto volume = SpectralVolumeData::fromMuVolume(volumeMu, /* std::move(volumeMu), if no longer required. */
database::attenuationModel(database::Composite::Bone_Cortical),
65.0f);
qInfo() << volume.max(); // output: 1.91896 (g/cm^3) <-- density representation
See also
muVolume()

◆ hasSpectralInformation()

bool CTL::SpectralVolumeData::hasSpectralInformation ( ) const

Returns true if this instance has full spectral information. This means the following conditions are fulfilled:

  • a data model describing the spectral dependency of the mass-attenuation coefficient for the material represented by this volume has been set
  • one of the following conditions is fulfilled:
    • numerical values in individual voxels represent material density (in g/cm^3), or
    • numerical values in individual voxels represent attenuation coefficiens (in 1/mm) and the corresponding reference energy is specified.

◆ isDensityVolume()

bool CTL::SpectralVolumeData::isDensityVolume ( ) const

Returns true if the data stored by this instance are density values (in g/cm^3).

auto volume = SpectralVolumeData(VoxelVolume<float>::cube(100, 1.0f, 0.02f),
database::attenuationModel(database::Composite::Water),
65.0f);
qInfo() << volume.isDensityVolume(); // output: false <-- constructor keeps attenuation domain
// This time, we use the factory method.
database::attenuationModel(database::Composite::Water),
65.0f);
qInfo() << volume2.isDensityVolume(); // output: true <-- fromMuVolume() transforms to density

◆ isMuVolume()

bool CTL::SpectralVolumeData::isMuVolume ( ) const

Returns true if the data stored by this instance are attenuation values (in 1/mm).

auto volume = SpectralVolumeData(VoxelVolume<float>::cube(100, 1.0f, 0.02f),
database::attenuationModel(database::Composite::Water),
65.0f);
qInfo() << volume.isDensityVolume(); // output: false <-- constructor keeps attenuation domain
// This time, we use the factory method.
database::attenuationModel(database::Composite::Water),
65.0f);
qInfo() << volume2.isDensityVolume(); // output: true <-- fromMuVolume() transforms to density

◆ massAttenuationCoeff()

float CTL::SpectralVolumeData::massAttenuationCoeff ( float  atEnergy) const

Returns the mass attenuation coefficient (in cm^2/g) of the material described by this instance for energy atEnergy (in keV).

Same as absorptionModel()->valueAt(atEnergy).

◆ materialName()

const QString & CTL::SpectralVolumeData::materialName ( ) const

Returns the name of the material described by this instance.

◆ meanMassAttenuationCoeff()

float CTL::SpectralVolumeData::meanMassAttenuationCoeff ( float  centerEnergy,
float  binWidth 
) const

Returns the mass attenuation coefficient of the material described by this instance averaged over the energy bin [centerEnergy - binWidth / 2, centerEnergy + binWidth / 2].

Same as absorptionModel()->meanValue(centerEnergy, binWidth).

◆ muVolume() [1/2]

std::unique_ptr< SpectralVolumeData > CTL::SpectralVolumeData::muVolume ( float  referenceEnergy) const

Returns the attenuation coefficient (with respect to the reference energy referenceEnergy) representation of this instance.

Note that this creates a copy of the data. The density values are transformed to attenuation coefficients with respect to the given reference energy referenceEnergy. In case this instance does already contain attenuation coefficient data, the values are re-referenced to referenceEnergy.

If this instance contains attenuation coefficient data and referenceEnergy == referenceEnergy(), it is discouraged to call this method because it would only create an unnecessary copy.

// Starting point: voxel data with attenuation coefficients in 1/mm
// here: make a cube phantom (100^3 voxels, size 1 mm^3) filled with value 0.02/mm (water)
VoxelVolume<float> muValues(100, 100, 100, 1.0f, 1.0f, 1.0f);
muValues.fill(0.02f); // 0.02 / mm is approx. attenuation of water for 60 keV
// Now: create spectral data by complementing \c muValues with the attenuation model of water
auto volume = SpectralVolumeData(muValues, /* std::move(muValues), if no longer required. */
database::attenuationModel(database::Composite::Water),
60.0f);
qInfo() << volume.max(); // output: 0.02 (our input, still in attenuation domain)
qInfo() << volume.referenceEnergy(); // output: 60
qInfo() << volume.referenceMassAttenuationCoeff(); // output: 0.2059
// auto otherVolume = volume.muVolume(60.0f); // meaningless copy, result identical to 'volume'
auto otherVolume = volume.muVolume(45.0f); // changes reference energy from 60 keV to 45 keV
qInfo() << otherVolume.max(); // output: 0.0240505 (rescaled to new reference energy)
qInfo() << otherVolume.referenceEnergy(); // output: 45
qInfo() << otherVolume.referenceMassAttenuationCoeff(); // output: 0.2476

◆ muVolume() [2/2]

std::unique_ptr< SpectralVolumeData > CTL::SpectralVolumeData::muVolume ( float  centerEnergy,
float  binWidth 
) const

Returns the attenuation coefficient with respect to an average mass attenuation coefficient in the energy interval [centerEnergy - binWidth / 2, centerEnergy + binWidth / 2].

Note that this creates a copy of the data.

See also
muVolume(float), meanMassAttenuationCoeff().

◆ referenceEnergy()

float CTL::SpectralVolumeData::referenceEnergy ( ) const

Returns the reference energy corresponding to the attenuation values managed by this instance. Does not contain meaningful information in case this instance manages density values.

◆ referenceMassAttenuationCoeff()

float CTL::SpectralVolumeData::referenceMassAttenuationCoeff ( ) const

Returns the reference mass attenuation coefficient corresponding to the attenuation values managed by this instance. Does not contain meaningful information in case this instance manages density values.

◆ replaceAbsorptionModel() [1/2]

void CTL::SpectralVolumeData::replaceAbsorptionModel ( AbstractIntegrableDataModel absorptionModel)

Replaces the absorption model in this instance by absorptionModel (must not be nullptr).

Note that this is only allowed for objects that already have spectral information available and throws an std::runtime_error otherwise.

◆ replaceAbsorptionModel() [2/2]

void CTL::SpectralVolumeData::replaceAbsorptionModel ( std::shared_ptr< AbstractIntegrableDataModel absorptionModel)

Replaces the absorption model in this instance by absorptionModel (must not be nullptr).

Note that this is only allowed for objects that already have spectral information available and throws an std::runtime_error otherwise.

◆ setDensity()

void CTL::SpectralVolumeData::setDensity ( VoxelVolume< float >  density)

Replaces the voxel data in this instance by the density values given by density (in g/cm^3).

Note that this is only allowed for objects that have a non-defaulted absorbtion model set and throws an std::runtime_error otherwise.

◆ setMaterialName()

void CTL::SpectralVolumeData::setMaterialName ( const QString &  name)

Sets the name of the material described by this instance to name.


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