Source code for openmc.source

from enum import Enum
from numbers import Real
from xml.etree import ElementTree as ET

import numpy as np
import h5py

import openmc.checkvalue as cv
from openmc.stats.multivariate import UnitSphere, Spatial
from openmc.stats.univariate import Univariate
from ._xml import get_text


[docs]class Source: """Distribution of phase space coordinates for source sites. Parameters ---------- space : openmc.stats.Spatial Spatial distribution of source sites angle : openmc.stats.UnitSphere Angular distribution of source sites energy : openmc.stats.Univariate Energy distribution of source sites filename : str Source file from which sites should be sampled library : str Path to a custom source library parameters : str Parameters to be provided to the custom source library .. versionadded:: 0.12 strength : float Strength of the source particle : {'neutron', 'photon'} Source particle type Attributes ---------- space : openmc.stats.Spatial or None Spatial distribution of source sites angle : openmc.stats.UnitSphere or None Angular distribution of source sites energy : openmc.stats.Univariate or None Energy distribution of source sites file : str or None Source file from which sites should be sampled library : str or None Path to a custom source library parameters : str Parameters to be provided to the custom source library strength : float Strength of the source particle : {'neutron', 'photon'} Source particle type """ def __init__(self, space=None, angle=None, energy=None, filename=None, library=None, parameters=None, strength=1.0, particle='neutron'): self._space = None self._angle = None self._energy = None self._file = None self._library = None self._parameters = None if space is not None: self.space = space if angle is not None: self.angle = angle if energy is not None: self.energy = energy if filename is not None: self.file = filename if library is not None: self.library = library if parameters is not None: self.parameters = parameters self.strength = strength self.particle = particle @property def file(self): return self._file @property def library(self): return self._library @property def parameters(self): return self._parameters @property def space(self): return self._space @property def angle(self): return self._angle @property def energy(self): return self._energy @property def strength(self): return self._strength @property def particle(self): return self._particle @file.setter def file(self, filename): cv.check_type('source file', filename, str) self._file = filename @library.setter def library(self, library_name): cv.check_type('library', library_name, str) self._library = library_name @parameters.setter def parameters(self, parameters_path): cv.check_type('parameters', parameters_path, str) self._parameters = parameters_path @space.setter def space(self, space): cv.check_type('spatial distribution', space, Spatial) self._space = space @angle.setter def angle(self, angle): cv.check_type('angular distribution', angle, UnitSphere) self._angle = angle @energy.setter def energy(self, energy): cv.check_type('energy distribution', energy, Univariate) self._energy = energy @strength.setter def strength(self, strength): cv.check_type('source strength', strength, Real) cv.check_greater_than('source strength', strength, 0.0, True) self._strength = strength @particle.setter def particle(self, particle): cv.check_value('source particle', particle, ['neutron', 'photon']) self._particle = particle
[docs] def to_xml_element(self): """Return XML representation of the source Returns ------- element : xml.etree.ElementTree.Element XML element containing source data """ element = ET.Element("source") element.set("strength", str(self.strength)) if self.particle != 'neutron': element.set("particle", self.particle) if self.file is not None: element.set("file", self.file) if self.library is not None: element.set("library", self.library) if self.parameters is not None: element.set("parameters", self.parameters) if self.space is not None: element.append(self.space.to_xml_element()) if self.angle is not None: element.append(self.angle.to_xml_element()) if self.energy is not None: element.append(self.energy.to_xml_element('energy')) return element
[docs] @classmethod def from_xml_element(cls, elem): """Generate source from an XML element Parameters ---------- elem : xml.etree.ElementTree.Element XML element Returns ------- openmc.Source Source generated from XML element """ source = cls() strength = get_text(elem, 'strength') if strength is not None: source.strength = float(strength) particle = get_text(elem, 'particle') if particle is not None: source.particle = particle filename = get_text(elem, 'file') if filename is not None: source.file = filename library = get_text(elem, 'library') if library is not None: source.library = library parameters = get_text(elem, 'parameters') if parameters is not None: source.parameters = parameters space = elem.find('space') if space is not None: source.space = Spatial.from_xml_element(space) angle = elem.find('angle') if angle is not None: source.angle = UnitSphere.from_xml_element(angle) energy = elem.find('energy') if energy is not None: source.energy = Univariate.from_xml_element(energy) return source
class ParticleType(Enum): NEUTRON = 0 PHOTON = 1 ELECTRON = 2 POSITRON = 3
[docs]class SourceParticle: """Source particle This class can be used to create source particles that can be written to a file and used by OpenMC Parameters ---------- r : iterable of float Position of particle in Cartesian coordinates u : iterable of float Directional cosines E : float Energy of particle in [eV] wgt : float Weight of the particle delayed_group : int Delayed group particle was created in (neutrons only) surf_id : int Surface ID where particle is at, if any. particle : ParticleType Type of the particle """ def __init__(self, r=(0., 0., 0.), u=(0., 0., 1.), E=1.0e6, wgt=1.0, delayed_group=0, surf_id=0, particle=ParticleType.NEUTRON): self.r = tuple(r) self.u = tuple(u) self.E = float(E) self.wgt = float(wgt) self.delayed_group = delayed_group self.surf_id = surf_id self.particle = particle
[docs] def to_tuple(self): """Return source particle attributes as a tuple Returns ------- tuple Source particle attributes """ return (self.r, self.u, self.E, self.wgt, self.delayed_group, self.surf_id, self.particle.value)
[docs]def write_source_file(source_particles, filename, **kwargs): """Write a source file using a collection of source particles Parameters ---------- source_particles : iterable of SourceParticle Source particles to write to file filename : str or path-like Path to source file to write **kwargs Keyword arguments to pass to :class:`h5py.File` See Also -------- openmc.SourceParticle """ # Create compound datatype for source particles pos_dtype = np.dtype([('x', '<f8'), ('y', '<f8'), ('z', '<f8')]) source_dtype = np.dtype([ ('r', pos_dtype), ('u', pos_dtype), ('E', '<f8'), ('wgt', '<f8'), ('delayed_group', '<i4'), ('surf_id', '<i4'), ('particle', '<i4'), ]) # Create array of source particles cv.check_iterable_type("source particles", source_particles, SourceParticle) arr = np.array([s.to_tuple() for s in source_particles], dtype=source_dtype) # Write array to file kwargs.setdefault('mode', 'w') with h5py.File(filename, **kwargs) as fh: fh.attrs['filetype'] = np.string_("source") fh.create_dataset('source_bank', data=arr, dtype=source_dtype)