An open source method of characteristics neutron transport code.
|
The Material class represents a unique material and its relevant nuclear data (i.e., multigroup cross-sections) for neutron transport. More...
#include "src/Material.h"
Public Member Functions | |
Material (int id=0, const char *name="") | |
Constructor sets the ID and unique ID for the Material. More... | |
virtual | ~Material () |
Destructor deletes all cross-section data structures from memory. | |
int | getId () const |
Return the Material's user-defined ID. More... | |
char * | getName () const |
Return the user-defined name of the Material. More... | |
double | getVolume () |
Return the aggregate volume/area of all instances of this Material. More... | |
int | getNumInstances () |
Return the number of instances of this Material in the Geometry. More... | |
int | getNumEnergyGroups () const |
Returns the number of energy groups for this Material's nuclear data. More... | |
FP_PRECISION * | getSigmaT () |
Return the array of the Material's total cross-sections. More... | |
FP_PRECISION | getMaxSigmaT () |
Return the maximum of the Material's total cross-sections. More... | |
FP_PRECISION * | getSigmaS () |
Return the array of the Material's scattering cross-section matrix. More... | |
FP_PRECISION * | getSigmaA () |
Return the array of the Material's absorption group cross-sections. More... | |
FP_PRECISION * | getSigmaF () |
Return the array of the Material's fission cross-sections. More... | |
FP_PRECISION * | getNuSigmaF () |
Return the array of the Material's fission cross-sections multiplied by nu . More... | |
FP_PRECISION * | getChi () |
Return the array of the Material's chi . More... | |
FP_PRECISION * | getFissionMatrix () |
Return the array of the Material's fission matrix. More... | |
FP_PRECISION | getSigmaTByGroup (int group) |
Get the Material's total cross section for some energy group. More... | |
FP_PRECISION | getSigmaSByGroup (int origin, int destination) |
Get the Material's scattering cross section for some energy group. More... | |
FP_PRECISION | getSigmaAByGroup (int group) |
Get the Material's absorption cross section for some energy group. More... | |
FP_PRECISION | getSigmaFByGroup (int group) |
Get the Material's fission cross section for some energy group. More... | |
FP_PRECISION | getNuSigmaFByGroup (int group) |
Get the Material's nu-fission cross section for some energy group. More... | |
FP_PRECISION | getChiByGroup (int group) |
Get the Material's fission spectrum for some energy group. More... | |
FP_PRECISION | getFissionMatrixByGroup (int origin, int destination) |
Get the Material's fission matrix for some energy group. More... | |
bool | isFissionable () |
Returns whether or not the Material contains a fissionable (non-zero) fission cross-section. More... | |
bool | isDataAligned () |
Returns true if the data is vector aligned, false otherwise (default). More... | |
int | getNumVectorGroups () |
Returns the rounded up number of energy groups to fill an integral number of vector lengths. More... | |
void | setName (const char *name) |
Sets the name of the Material. More... | |
void | setVolume (double volume) |
Set the volume/area of the Material. More... | |
void | incrementVolume (double volume) |
Increment the volume/area of the Material by some amount. More... | |
void | setNumInstances (int num_instances) |
Set the number of instances of this Material. More... | |
void | incrementNumInstances () |
Increment the number of instances of this Material. More... | |
void | setNumEnergyGroups (const int num_groups) |
Set the number of energy groups for this Material. More... | |
void | setSigmaT (double *xs, int num_groups) |
Set the Material's array of total cross-sections. More... | |
void | setSigmaS (double *xs, int num_groups) |
Set the Material's 2D array of scattering cross-sections. More... | |
void | setSigmaF (double *xs, int num_groups) |
Set the Material's array of fission cross-sections. More... | |
void | setNuSigmaF (double *xs, int num_groups) |
Set the Material's array of fission cross-sections multiplied by . More... | |
void | setChi (double *xs, int num_groups) |
Set the Material's array of chi values. More... | |
void | setSigmaTByGroup (double xs, int group) |
Set the Material's total cross-section for some energy group. More... | |
void | setSigmaFByGroup (double xs, int group) |
Set the Material's fission cross-section for some energy group. More... | |
void | setNuSigmaFByGroup (double xs, int group) |
Set the Material's fission cross-section multiplied by for some energy group. More... | |
void | setSigmaSByGroup (double xs, int origin, int destination) |
Set the Material's scattering cross-section for some energy group. More... | |
void | setChiByGroup (double xs, int group) |
Set the Material's chi value for some energy group. More... | |
void | setSigmaAByGroup (double xs, int group) |
Set the Material's absorption cross-section for some energy group. More... | |
void | buildFissionMatrix () |
Builds the fission matrix from chi and the fission cross-section. More... | |
void | transposeProductionMatrices () |
Transposes the scattering and fission matrices. More... | |
void | alignData () |
Reallocates the Material's cross-section data structures along word-aligned boundaries. More... | |
Material * | clone () |
Create a duplicate of the Material. More... | |
std::string | toString () |
Converts this Material's attributes to a character array representation. More... | |
void | printString () |
Prints a string representation of all of the Material's attributes to the console. | |
The Material class represents a unique material and its relevant nuclear data (i.e., multigroup cross-sections) for neutron transport.
Material::Material | ( | int | id = 0 , |
const char * | name = "" |
||
) |
void Material::alignData | ( | ) |
Reallocates the Material's cross-section data structures along word-aligned boundaries.
This method is used to assist with SIMD auto-vectorization of the MOC routines in the Solver classes. Rather than using the assigned number of energy groups, this method adds "dummy" energy groups such that the total number of groups is some multiple of VEC_LENGTH (typically 4, 8, or 16). As a result, the SIMD-vectorized Solver subclasses can break up loops over energy groups in such a way to "expose" the SIMD nature of the algorithm.
void Material::buildFissionMatrix | ( | ) |
Builds the fission matrix from chi and the fission cross-section.
The fission matrix is constructed as the outer product of the chi and fission cross-section vectors. This routine is intended for internal use and is called by the Solver at runtime.
FP_PRECISION * Material::getChi | ( | ) |
FP_PRECISION Material::getChiByGroup | ( | int | group | ) |
Get the Material's fission spectrum for some energy group.
group | the energy group |
FP_PRECISION * Material::getFissionMatrix | ( | ) |
FP_PRECISION Material::getFissionMatrixByGroup | ( | int | origin, |
int | destination | ||
) |
Get the Material's fission matrix for some energy group.
origin | the incoming energy group |
destination | the outgoing energy group |
int Material::getId | ( | ) | const |
FP_PRECISION Material::getMaxSigmaT | ( | ) |
char * Material::getName | ( | ) | const |
int Material::getNumEnergyGroups | ( | ) | const |
Returns the number of energy groups for this Material's nuclear data.
int Material::getNumInstances | ( | ) |
int Material::getNumVectorGroups | ( | ) |
Returns the rounded up number of energy groups to fill an integral number of vector lengths.
FP_PRECISION * Material::getNuSigmaF | ( | ) |
FP_PRECISION Material::getNuSigmaFByGroup | ( | int | group | ) |
Get the Material's nu-fission cross section for some energy group.
group | the energy group |
FP_PRECISION * Material::getSigmaA | ( | ) |
FP_PRECISION Material::getSigmaAByGroup | ( | int | group | ) |
Get the Material's absorption cross section for some energy group.
group | the energy group |
FP_PRECISION * Material::getSigmaF | ( | ) |
FP_PRECISION Material::getSigmaFByGroup | ( | int | group | ) |
Get the Material's fission cross section for some energy group.
group | the energy group |
FP_PRECISION * Material::getSigmaS | ( | ) |
FP_PRECISION Material::getSigmaSByGroup | ( | int | origin, |
int | destination | ||
) |
Get the Material's scattering cross section for some energy group.
origin | the incoming energy group |
destination | the outgoing energy group |
FP_PRECISION * Material::getSigmaT | ( | ) |
FP_PRECISION Material::getSigmaTByGroup | ( | int | group | ) |
Get the Material's total cross section for some energy group.
group | the energy group |
double Material::getVolume | ( | ) |
void Material::incrementNumInstances | ( | ) |
Increment the number of instances of this Material.
This routine is called by the TrackGenerator during track generation and segmentation.
void Material::incrementVolume | ( | double | volume | ) |
Increment the volume/area of the Material by some amount.
This routine is called by the TrackGenerator during track generation and segmentation.
volume | the amount to increment the current volume by |
bool Material::isDataAligned | ( | ) |
Returns true if the data is vector aligned, false otherwise (default).
bool Material::isFissionable | ( | ) |
Returns whether or not the Material contains a fissionable (non-zero) fission cross-section.
void Material::setChi | ( | double * | xs, |
int | num_groups | ||
) |
Set the Material's array of chi values.
This method is a helper function to allow OpenMOC users to assign the Material's nuclear data in Python. A user must initialize a NumPy array of the correct size (e.g., a float64 array the length of the number of energy groups) as input to this function. This function then fills the NumPy array with the data values for the Material's chi distribution. An example of how this function might be called in Python is as follows:
xs | the array of chi values |
num_groups | the number of energy groups |
void Material::setChiByGroup | ( | double | xs, |
int | group | ||
) |
Set the Material's chi value for some energy group.
xs | the chi value ( ) |
group | the energy group |
void Material::setName | ( | const char * | name | ) |
void Material::setNumEnergyGroups | ( | const int | num_groups | ) |
Set the number of energy groups for this Material.
num_groups | the number of energy groups. |
void Material::setNumInstances | ( | int | num_instances | ) |
void Material::setNuSigmaF | ( | double * | xs, |
int | num_groups | ||
) |
Set the Material's array of fission cross-sections multiplied by .
xs | the array of fission cross-sections multiplied by nu |
num_groups | the number of energy groups |
void Material::setNuSigmaFByGroup | ( | double | xs, |
int | group | ||
) |
Set the Material's fission cross-section multiplied by for some energy group.
This method is a helper function to allow OpenMOC users to assign the Material's nuclear data in Python. A user must initialize a NumPy array of the correct size (e.g., a float64 array the length of the number of energy groups) as input to this function. This function then fills the NumPy array with the data values for the Material's nu*fission cross-sections. An example of how this function might be called in Python is as follows:
xs | the fission cross-section multiplied by nu |
group | the energy group |
void Material::setSigmaAByGroup | ( | double | xs, |
int | group | ||
) |
Set the Material's absorption cross-section for some energy group.
The absorption cross section is computed by OpenMOC when needed, but the user can set it manually, to replace absorption by another cross section as absorption is not needed when determining Keff from fission rates.
xs | the absorption cross-section |
group | the energy group |
void Material::setSigmaF | ( | double * | xs, |
int | num_groups | ||
) |
Set the Material's array of fission cross-sections.
This method is a helper function to allow OpenMOC users to assign the Material's nuclear data in Python. A user must initialize a NumPy array of the correct size (e.g., a float64 array the length of the number of energy groups) as input to this function. This function then fills the NumPy array with the data values for the Material's fission cross-sections. An example of how this function might be called in Python is as follows:
xs | the array of fission cross-sections |
num_groups | the number of energy groups |
void Material::setSigmaFByGroup | ( | double | xs, |
int | group | ||
) |
Set the Material's fission cross-section for some energy group.
xs | the fission cross-section |
group | the energy group |
void Material::setSigmaS | ( | double * | xs, |
int | num_groups_squared | ||
) |
Set the Material's 2D array of scattering cross-sections.
The array should be passed to OpenMOC as a 1D array in column-major order. This assumes the standard convention, where column index is the origin group and the row index is the destination group. That is, the array should be ordered as follows: 1 -> 1 1 -> 2 1 -> 3 ... 2 -> 1 2 -> 2 ...
Note that if the scattering matrix is defined in NumPy by the standard convention, "flat" will put the matrix into row major order. Thus, one should transpose the matrix before flattening.
For cache efficiency, the transpose of the input is actually stored in OpenMOC.
This method is a helper function to allow OpenMOC users to assign the Material's nuclear data in Python. A user must initialize a NumPy array of the correct size (e.g., a float64 array the length of the square of the number of energy groups) as input to this function. This function then fills the NumPy array with the data values for the Material's scattering cross-sections. An example of how this function might be called in Python is as follows:
xs | the array of scattering cross-sections |
num_groups_squared | the number of energy groups squared |
void Material::setSigmaSByGroup | ( | double | xs, |
int | origin, | ||
int | destination | ||
) |
Set the Material's scattering cross-section for some energy group.
xs | the scattering cross-section |
origin | the column index in the scattering matrix |
destination | the row index in the scattering matrix |
void Material::setSigmaT | ( | double * | xs, |
int | num_groups | ||
) |
Set the Material's array of total cross-sections.
This method is a helper function to allow OpenMOC users to assign the Material's nuclear data in Python. A user must initialize a NumPy array of the correct size (e.g., a float64 array the length of the number of energy groups) as input to this function. This function then fills the NumPy array with the data values for the Material's total cross-sections. An example of how this function might be called in Python is as follows:
NOTE: This routine will override any zero-valued cross-sections (e.g., in void or gap regions) with a minimum value of 1E-10 to avoid numerical issues in the MOC solver.
xs | the array of total cross-sections |
num_groups | the number of energy groups |
void Material::setSigmaTByGroup | ( | double | xs, |
int | group | ||
) |
Set the Material's total cross-section for some energy group.
xs | the total cross-section |
group | the energy group |
void Material::setVolume | ( | double | volume | ) |
std::string Material::toString | ( | ) |
Converts this Material's attributes to a character array representation.
The character array returned includes the user-defined ID, and each of the absorption, total, fission, nu multiplied by fission and scattering cross-sections and chi for all energy groups.
void Material::transposeProductionMatrices | ( | ) |
Transposes the scattering and fission matrices.
This routine is used by the Solver when performing adjoint flux caclulations.