openmc.model.Model
- class openmc.model.Model(geometry: Geometry | None = None, materials: Materials | None = None, settings: Settings | None = None, tallies: Tallies | None = None, plots: Plots | None = None)[source]
Model container.
This class can be used to store instances of
openmc.Geometry,openmc.Materials,openmc.Settings,openmc.Tallies, andopenmc.Plots, thus making a complete model. TheModel.export_to_xml()method will export XML files for all attributes that have been set. If theModel.materialsattribute is not set, it will attempt to create amaterials.xmlfile based on all materials appearing in the geometry.Changed in version 0.13.0: The model information can now be loaded in to OpenMC directly via openmc.lib
- Parameters:
geometry (openmc.Geometry, optional) – Geometry information
materials (openmc.Materials, optional) – Materials information
settings (openmc.Settings, optional) – Settings information
tallies (openmc.Tallies, optional) – Tallies information
plots (openmc.Plots, optional) – Plot information
- Variables:
geometry (openmc.Geometry) – Geometry information
materials (openmc.Materials) – Materials information
settings (openmc.Settings) – Settings information
tallies (openmc.Tallies) – Tallies information
plots (openmc.Plots) – Plot information
- add_kinetics_parameters_tallies(num_groups: int | None = None)[source]
Add tallies for calculating kinetics parameters using the IFP method.
This method adds tallies to the model for calculating two kinetics parameters, the generation time and the effective delayed neutron fraction (beta effective). After a model is run, these parameters can be determined through the
openmc.StatePoint.ifp_results()method.- Parameters:
num_groups (int, optional) – Number of precursor groups to filter the delayed neutron fraction. If None, only the total effective delayed neutron fraction is tallied.
- apply_tally_results(statepoint: str | PathLike | StatePoint)[source]
Apply results from a statepoint to tally objects on the Model
- Parameters:
statepoint (PathLike or openmc.StatePoint) – Statepoint file used to update tally results
- calculate_volumes(threads: int | None = None, output: bool = True, cwd: str | PathLike = '.', openmc_exec: str | PathLike = 'openmc', mpi_args: list[str] | None = None, apply_volumes: bool = True, export_model_xml: bool = True, **export_kwargs)[source]
Runs an OpenMC stochastic volume calculation and, if requested, applies volumes to the model
Added in version 0.13.0.
- Parameters:
threads (int, optional) – Number of OpenMP threads. If OpenMC is compiled with OpenMP threading enabled, the default is implementation-dependent but is usually equal to the number of hardware threads available (or a value set by the
OMP_NUM_THREADSenvironment variable). This currenty only applies to the case when not using the C API.output (bool, optional) – Capture OpenMC output from standard out
openmc_exec (str, optional) – Path to OpenMC executable. Defaults to ‘openmc’. This only applies to the case when not using the C API.
mpi_args (list of str, optional) – MPI execute command and any additional MPI arguments to pass, e.g. [‘mpiexec’, ‘-n’, ‘8’]. This only applies to the case when not using the C API.
cwd (str, optional) – Path to working directory to run in. Defaults to the current working directory.
apply_volumes (bool, optional) – Whether apply the volume calculation results from this calculation to the model. Defaults to applying the volumes.
export_model_xml (bool, optional) – Exports a single model.xml file rather than separate files. Defaults to True.
**export_kwargs – Keyword arguments passed to either
Model.export_to_model_xml()orModel.export_to_xml().
- convert_to_multigroup(method: str = 'material_wise', groups: str = 'CASMO-2', nparticles: int = 2000, overwrite_mgxs_library: bool = False, mgxs_path: str | PathLike = 'mgxs.h5', correction: str | None = None)[source]
Convert all materials from continuous energy to multigroup.
If no MGXS data library file is found, generate one using one or more continuous energy Monte Carlo simulations.
- Parameters:
method ({"material_wise", "stochastic_slab", "infinite_medium"}, optional) – Method to generate the MGXS.
groups (openmc.mgxs.EnergyGroups or str, optional) – Energy group structure for the MGXS or the name of the group structure (based on keys from openmc.mgxs.GROUP_STRUCTURES).
mgxs_path (str, optional) – Filename of the mgxs.h5 library file.
correction (str, optional) – Transport correction to apply to the MGXS. Options are None and “P0”.
- convert_to_random_ray()[source]
Convert a multigroup model to use random ray.
This method determines values for the needed settings and adds them to the settings.random_ray dictionary so as to enable random ray mode. The settings that are populated are:
‘ray_source’ (openmc.IndependentSource): Where random ray starting points are sampled from.
‘distance_inactive’ (float): The “dead zone” distance at the beginning of the ray.
‘distance_active’ (float): The “active” distance of the ray
‘particles’ (int): Number of rays to simulate
The method will determine reasonable defaults for each of the above variables based on analysis of the model’s geometry. The function will have no effect if the random ray dictionary is already defined in the model settings.
- deplete(method: str = 'cecm', final_step: bool = True, operator_kwargs: dict | None = None, directory: str | PathLike = '.', output: bool = True, **integrator_kwargs)[source]
Deplete model using specified timesteps/power
Changed in version 0.13.0: The final_step, operator_kwargs, directory, and output arguments were added.
- Parameters:
timesteps (iterable of float or iterable of tuple) – Array of timesteps. Note that values are not cumulative. The units are specified by the timestep_units argument when timesteps is an iterable of float. Alternatively, units can be specified for each step by passing an iterable of (value, unit) tuples.
method (str) – Integration method used for depletion (e.g., ‘cecm’, ‘predictor’). Defaults to ‘cecm’.
final_step (bool, optional) – Indicate whether or not a transport solve should be run at the end of the last timestep. Defaults to running this transport solve.
operator_kwargs (dict) – Keyword arguments passed to the depletion operator initializer (e.g.,
openmc.deplete.Operator())directory (PathLike, optional) – Directory to write XML files to. If it doesn’t exist already, it will be created. Defaults to the current working directory
output (bool) – Capture OpenMC output from standard out
integrator_kwargs (dict) – Remaining keyword arguments passed to the depletion integrator (e.g.,
openmc.deplete.CECMIntegrator).
- differentiate_depletable_mats(diff_volume_method: str = None)[source]
Assign distribmats for each depletable material
Added in version 0.14.0.
Changed in version 0.15.1: diff_volume_method default is None, do not set volumes on the new material ovjects. Is now a convenience method for differentiate_mats(diff_volume_method, depletable_only=True)
- Parameters:
diff_volume_method (str) – Specifies how the volumes of the new materials should be found. - None: Do not assign volumes to the new materials (Default) - ‘divide equally’: Divide the original material volume equally between the new materials - ‘match cell’: Set the volume of the material to the volume of the cell they fill
- differentiate_mats(diff_volume_method: str = None, depletable_only: bool = True)[source]
Assign distribmats for each material
Added in version 0.15.1.
- Parameters:
diff_volume_method (str) – Specifies how the volumes of the new materials should be found. - None: Do not assign volumes to the new materials (Default) - ‘divide equally’: Divide the original material volume equally between the new materials - ‘match cell’: Set the volume of the material to the volume of the cell they fill
depletable_only (bool) – Default is True, only depletable materials will be differentiated. If False, all materials will be differentiated.
- export_to_model_xml(path: str | PathLike = 'model.xml', remove_surfs: bool = False, nuclides_to_ignore: Iterable[str] | None = None)[source]
Export model to a single XML file.
Added in version 0.13.3.
- Parameters:
- export_to_xml(directory: str | PathLike = '.', remove_surfs: bool = False, nuclides_to_ignore: Iterable[str] | None = None)[source]
Export model to separate XML files.
- Parameters:
directory (PathLike) – Directory to write XML files to. If it doesn’t exist already, it will be created.
remove_surfs (bool) –
Whether or not to remove redundant surfaces from the geometry when exporting.
Added in version 0.13.1.
nuclides_to_ignore (list of str) – Nuclides to ignore when exporting to XML.
- finalize_lib()[source]
Finalize simulation and free memory allocated for the C API
Added in version 0.13.0.
- classmethod from_model_xml(path: str | PathLike = 'model.xml') Model[source]
Create model from single XML file
Added in version 0.13.3.
- Parameters:
path (PathLike) – Path to model.xml file
- classmethod from_xml(geometry: str | PathLike = 'geometry.xml', materials: str | PathLike = 'materials.xml', settings: str | PathLike = 'settings.xml', tallies: str | PathLike = 'tallies.xml', plots: str | PathLike = 'plots.xml') Model[source]
Create model from existing XML files
- Parameters:
geometry (PathLike) – Path to geometry.xml file
materials (PathLike) – Path to materials.xml file
settings (PathLike) – Path to settings.xml file
tallies (PathLike) –
Path to tallies.xml file
Added in version 0.13.0.
plots (PathLike) –
Path to plots.xml file
Added in version 0.13.0.
- Returns:
Model created from XML files
- Return type:
- id_map(origin: Sequence[float] | None = None, width: Sequence[float] | None = None, pixels: int | Sequence[int] = 40000, basis: str = 'xy', **init_kwargs) ndarray[source]
Generate an ID map for domains based on the plot parameters
If the model is not yet initialized, it will be initialized with openmc.lib. If the model is initialized, the model will remain initialized after this method call exits.
Added in version 0.15.3.
- Parameters:
origin (Sequence[float], optional) – Origin of the plot. If unspecified, this argument defaults to the center of the bounding box if the bounding box does not contain inf values for the provided basis, otherwise (0.0, 0.0, 0.0).
width (Sequence[float], optional) – Width of the plot. If unspecified, this argument defaults to the width of the bounding box if the bounding box does not contain inf values for the provided basis, otherwise (10.0, 10.0).
pixels (int | Sequence[int], optional) – If an iterable of ints is provided then this directly sets the number of pixels to use in each basis direction. If a single int is provided then this sets the total number of pixels in the plot and the number of pixels in each basis direction is calculated from this total and the image aspect ratio based on the width argument.
basis ({'xy', 'yz', 'xz'}, optional) – Basis of the plot.
**init_kwargs – Keyword arguments passed to
Model.init_lib().
- Returns:
id_map – A NumPy array with shape (vertical pixels, horizontal pixels, 3) of OpenMC property IDs with dtype int32. The last dimension of the array contains cell IDs, cell instances, and material IDs (in that order).
- Return type:
- import_properties(filename: str | PathLike)[source]
Import physical properties
Changed in version 0.13.0: This method now updates values as loaded in memory with the C API
- Parameters:
filename (PathLike) – Path to properties HDF5 file
See also
- init_lib(threads: int | None = None, geometry_debug: bool = False, restart_file: str | PathLike | None = None, tracks: bool = False, output: bool = True, event_based: bool | None = None, intracomm=None, directory: str | PathLike | None = None)[source]
Initializes the model in memory via the C API
Added in version 0.13.0.
- Parameters:
threads (int, optional) – Number of OpenMP threads. If OpenMC is compiled with OpenMP threading enabled, the default is implementation-dependent but is usually equal to the number of hardware threads available (or a value set by the
OMP_NUM_THREADSenvironment variable).geometry_debug (bool, optional) – Turn on geometry debugging during simulation. Defaults to False.
restart_file (PathLike, optional) – Path to restart file to use
tracks (bool, optional) – Enables the writing of particles tracks. The number of particle tracks written to tracks.h5 is limited to 1000 unless Settings.max_tracks is set. Defaults to False.
output (bool) – Capture OpenMC output from standard out
event_based (None or bool, optional) – Turns on event-based parallelism if True. If None, the value in the Settings will be used.
intracomm (mpi4py.MPI.Intracomm or None, optional) – MPI intracommunicator
directory (PathLike or None, optional) – Directory to write XML files to. Defaults to None.
- keff_search(func: ModelModifier, x0: float, x1: float, target: float = 1.0, k_tol: float = 0.0001, sigma_final: float = 0.0003, p: float = 0.5, q: float = 0.95, memory: int = 4, x_min: float | None = None, x_max: float | None = None, b0: int | None = None, b_min: int = 20, b_max: int | None = None, maxiter: int = 50, output: bool = False, func_kwargs: dict[str, Any] | None = None, run_kwargs: dict[str, Any] | None = None) SearchResult[source]
Perform a keff search on a model parametrized by a single variable.
This method uses the GRsecant method described in a paper by Price and Roskoff. The GRsecant method is a modification of the secant method that accounts for uncertainties in the function evaluations. The method uses a weighted linear fit of the most recent function evaluations to predict the next point to evaluate. It also adaptively changes the number of batches to meet the target uncertainty value at each iteration.
The target uncertainty for iteration \(n+1\) is determined by the following equation (following Eq. (8) in the paper):
\[\sigma_{i+1} = q \sigma_\text{final} \left ( \frac{ \min \left \{ \left\lvert k_i - k_\text{target} \right\rvert : k=0,1,\dots,n \right \} }{k_\text{tol}} \right )^p \]where \(q\) is a multiplicative factor less than 1, given as the
sigma_factorparameter below.- Parameters:
func (ModelModifier) – Function that takes the parameter to be searched and makes a modification to the model.
x0 (float) – First guess for the parameter passed to func
x1 (float) – Second guess for the parameter passed to func
target (float, optional) – keff value to search for
k_tol (float, optional) – Stopping criterion on the function value; the absolute value must be within
k_tolof zero to be accepted.sigma_final (float, optional) – Maximum accepted k-effective uncertainty for the stopping criterion.
p (float, optional) – Exponent used in the stopping criterion.
q (float, optional) – Multiplicative factor used in the stopping criterion.
memory (int, optional) – Number of most-recent points used in the weighted linear fit of
f(x) = a + b xto predict the next point.x_min (float, optional) – Minimum allowed value for the parameter
x.x_max (float, optional) – Maximum allowed value for the parameter
x.b0 (int, optional) – Number of active batches to use for the initial function evaluations. If None, uses the model’s current setting.
b_min (int, optional) – Minimum number of active batches to use in a function evaluation.
b_max (int, optional) – Maximum number of active batches to use in a function evaluation.
maxiter (int, optional) – Maximum number of iterations to perform.
output (bool, optional) – Whether or not to display output showing iteration progress.
func_kwargs (dict, optional) – Keyword-based arguments to pass to the func function.
run_kwargs (dict, optional) – Keyword arguments to pass to
openmc.Model.run()oropenmc.lib.run().
- Returns:
Result object containing the estimated root (parameter value) and evaluation history (parameters, means, standard deviations, and batches), plus convergence status and termination reason.
- Return type:
SearchResult
- plot(origin: Sequence[float] | None = None, width: Sequence[float] | None = None, pixels: int | Sequence[int] = 40000, basis: str = 'xy', color_by: str = 'cell', colors: dict | None = None, seed: int | None = None, openmc_exec: str | PathLike = 'openmc', axes=None, legend: bool = False, axis_units: str = 'cm', outline: bool | str = False, show_overlaps: bool = False, overlap_color: Sequence[int] | str | None = None, n_samples: int | None = None, plane_tolerance: float = 1.0, legend_kwargs: dict | None = None, source_kwargs: dict | None = None, contour_kwargs: dict | None = None, **kwargs)[source]
Display a slice plot of the model.
Added in version 0.15.1.
- Parameters:
origin (iterable of float) – Coordinates at the origin of the plot. If left as None, the center of the bounding box will be used to attempt to ascertain the origin with infinite values being replaced by 0.
width (iterable of float) – Width of the plot in each basis direction. If left as none then the width of the bounding box will be used to attempt to ascertain the plot width. Defaults to (10, 10) if the bounding box contains inf values.
pixels (Iterable of int or int) – If an iterable of ints is provided then this directly sets the number of pixels to use in each basis direction. If a single int is provided then this sets the total number of pixels in the plot and the number of pixels in each basis direction is calculated from this total and the image aspect ratio based on the width argument.
basis ({'xy', 'xz', 'yz'}) – The basis directions for the plot
color_by ({'cell', 'material'}) – Indicate whether the plot should be colored by cell or by material
colors (dict) –
Assigns colors to specific materials or cells. Keys are instances of
CellorMaterialand values are RGB 3-tuples, RGBA 4-tuples, or strings indicating SVG color names. Red, green, blue, and alpha should all be floats in the range [0.0, 1.0], for example:# Make water blue water = openmc.Cell(fill=h2o) universe.plot(..., colors={water: (0., 0., 1.))
seed (int) – Seed for the random number generator
openmc_exec (str) – Path to OpenMC executable.
axes (matplotlib.Axes) –
Axes to draw to
Added in version 0.13.1.
legend (bool) –
Whether a legend showing material or cell names should be drawn
Added in version 0.14.0.
axis_units ({'km', 'm', 'cm', 'mm'}) –
Units used on the plot axis
Added in version 0.14.0.
Whether outlines between color boundaries should be drawn. If set to ‘only’, only outlines will be drawn.
Added in version 0.14.0.
show_overlaps (bool) – Indicate whether or not overlapping regions are shown. Default is False.
overlap_color (Iterable of int or str) – Color to apply to overlapping regions. Default is red.
n_samples (int, optional) – The number of source particles to sample and add to plot. Defaults to None which doesn’t plot any particles on the plot.
plane_tolerance (float) – When plotting a plane the source locations within the plane +/- the plane_tolerance will be included and those outside of the plane_tolerance will not be shown
legend_kwargs (dict) –
Keyword arguments passed to
matplotlib.pyplot.legend().Added in version 0.14.0.
source_kwargs (dict, optional) – Keyword arguments passed to
matplotlib.pyplot.scatter().contour_kwargs (dict, optional) – Keyword arguments passed to
matplotlib.pyplot.contour().**kwargs – Keyword arguments passed to
matplotlib.pyplot.imshow().
- Returns:
Axes containing resulting image
- Return type:
- plot_geometry(output: bool = True, cwd: str | PathLike = '.', openmc_exec: str | PathLike = 'openmc', export_model_xml: bool = True, **export_kwargs)[source]
Creates plot images as specified by the Model.plots attribute
Added in version 0.13.0.
- Parameters:
output (bool, optional) – Capture OpenMC output from standard out
cwd (PathLike, optional) – Path to working directory to run in. Defaults to the current working directory.
openmc_exec (PathLike, optional) – Path to OpenMC executable. Defaults to ‘openmc’. This only applies to the case when not using the C API.
export_model_xml (bool, optional) – Exports a single model.xml file rather than separate files. Defaults to True.
**export_kwargs – Keyword arguments passed to either
Model.export_to_model_xml()orModel.export_to_xml().
- rotate_cells(names_or_ids: Iterable[str] | Iterable[int], vector: Iterable[float])[source]
Rotate the identified cell(s) by the specified rotation vector. The rotation is only applied to cells filled with a universe.
Note
If applying this change to a name that is not unique, then the change will be applied to all objects of that name.
Added in version 0.13.0.
- Parameters:
names_or_ids (Iterable of str or int) – The cell names (if str) or id (if int) that are to be translated or rotated. This parameter can include a mix of names and ids.
vector (Iterable of float) – The rotation vector of length 3 to apply. This array specifies the angles in degrees about the x, y, and z axes, respectively.
- run(particles: int | None = None, threads: int | None = None, geometry_debug: bool = False, restart_file: str | PathLike | None = None, tracks: bool = False, output: bool = True, cwd: str | PathLike = '.', openmc_exec: str | PathLike = 'openmc', mpi_args: Iterable[str] = None, event_based: bool | None = None, export_model_xml: bool = True, apply_tally_results: bool = False, **export_kwargs) Path[source]
Run OpenMC
If the C API has been initialized, then the C API is used, otherwise, this method creates the XML files and runs OpenMC via a system call. In both cases this method returns the path to the last statepoint file generated.
Changed in version 0.12: Instead of returning the final k-effective value, this function now returns the path to the final statepoint written.
Changed in version 0.13.0: This method can utilize the C API for execution
- Parameters:
particles (int, optional) – Number of particles to simulate per generation.
threads (int, optional) – Number of OpenMP threads. If OpenMC is compiled with OpenMP threading enabled, the default is implementation-dependent but is usually equal to the number of hardware threads available (or a value set by the
OMP_NUM_THREADSenvironment variable).geometry_debug (bool, optional) – Turn on geometry debugging during simulation. Defaults to False.
restart_file (str or PathLike) – Path to restart file to use
tracks (bool, optional) – Enables the writing of particles tracks. The number of particle tracks written to tracks.h5 is limited to 1000 unless Settings.max_tracks is set. Defaults to False.
output (bool, optional) – Capture OpenMC output from standard out
cwd (PathLike, optional) – Path to working directory to run in. Defaults to the current working directory.
openmc_exec (str, optional) – Path to OpenMC executable. Defaults to ‘openmc’.
mpi_args (list of str, optional) – MPI execute command and any additional MPI arguments to pass, e.g. [‘mpiexec’, ‘-n’, ‘8’].
event_based (None or bool, optional) – Turns on event-based parallelism if True. If None, the value in the Settings will be used.
export_model_xml (bool, optional) –
Exports a single model.xml file rather than separate files. Defaults to True.
Added in version 0.13.3.
apply_tally_results (bool) –
Whether to apply results of the final statepoint file to the model’s tally objects.
Added in version 0.15.1.
**export_kwargs – Keyword arguments passed to either
Model.export_to_model_xml()orModel.export_to_xml().
- Returns:
Path to the last statepoint written by this run (None if no statepoint was written)
- Return type:
Path
- sample_external_source(n_samples: int = 1000, prn_seed: int | None = None, **init_kwargs) ParticleList[source]
Sample external source and return source particles.
Added in version 0.15.1.
- Parameters:
n_samples (int) – Number of samples
prn_seed (int) – Pseudorandom number generator (PRNG) seed; if None, one will be generated randomly.
**init_kwargs – Keyword arguments passed to
openmc.lib.init()
- Returns:
List of samples source particles
- Return type:
- sync_dagmc_universes()[source]
Synchronize all DAGMC universes in the current geometry.
This method iterates over all DAGMC universes in the geometry and synchronizes their cells with the current material assignments. Requires that the model has been initialized via
Model.init_lib().Added in version 0.15.1.
- translate_cells(names_or_ids: Iterable[str] | Iterable[int], vector: Iterable[float])[source]
Translate the identified cell(s) by the specified translation vector. The translation is only applied to cells filled with a universe.
Note
If applying this change to a name that is not unique, then the change will be applied to all objects of that name.
Added in version 0.13.0.
- Parameters:
names_or_ids (Iterable of str or int) – The cell names (if str) or id (if int) that are to be translated or rotated. This parameter can include a mix of names and ids.
vector (Iterable of float) – The translation vector of length 3 to apply. This array specifies the x, y, and z dimensions of the translation.
- update_cell_temperatures(names_or_ids: Iterable[str] | Iterable[int], temperature: float)[source]
Update the temperature of a set of cells to the given value
Note
If applying this change to a name that is not unique, then the change will be applied to all objects of that name.
Added in version 0.13.0.
- update_densities(names_or_ids: Iterable[str] | Iterable[int], density: float, density_units: str = 'atom/b-cm')[source]
Update the density of a given set of materials to a new value
Note
If applying this change to a name that is not unique, then the change will be applied to all objects of that name.
Added in version 0.13.0.
- Parameters:
names_or_ids (Iterable of str or int) – The material names (if str) or id (if int) that are to be updated. This parameter can include a mix of names and ids.
density (float) – The density to apply in the units specified by density_units
density_units ({'atom/b-cm', 'g/cm3'}, optional) – Units for density. Defaults to ‘atom/b-cm’