sandy.core.endf6 module
Created on Wed Dec 4 14:50:33 2019
@author: lfiorito
- class sandy.core.endf6.Endf6(data, file=None)
Bases:
_FormattedFile
Container for ENDF-6 file text grouped by MAT, MF and MT numbers.
- Attributes:
- data
- is_empty
keys
List of keys (MAT, MF, MT) used to identify each tape section.
kind
Kind of ENDF-6 formatted file (‘endf6’, ‘pendf’, ‘gendf’, ‘errorr’) .
- mat
- mf
- mt
Methods
get_ace
([temperature, err, ...])- Parameters:
get_pendf
([temperature, err, ...])- Parameters:
get_errorr
([temperature, err, ...])- Parameters:
get_id
([method])Extract ID for a given MAT for a ENDF-6 file.
get_nsub
()Determine ENDF-6 sub-library type by reading flag "NSUB" of first MAT in file.
Extract MAT, MF and MT combinations avaialbel in the file and report it in tabulated format.
read_section
(mat, mf, mt[, raise_error])Parse MAT/MF/MT section.
update_intro
(**kwargs)Method to update MF1/MT451 of each MAT based on the file content (concistency is enforced) and user-given keyword arguments.
- apply_perturbations(smps, processes=1, njoy_kws={}, **kwargs)
Apply perturbations to the data contained in ENDF6 file. At the moment only the procedure for cross sections is implemented. Options are included to directly convert perturbed pendf to ace and write data on files.
- Parameters:
- smpssamples obtained taking the relative covariances from the
- evaluated files and a unit vector as mean.
- processesnumber of processes used to complete the task.
Creation of perturbed pendf files and conversion to ACE format is done in parallel if processes>1. The default is 1.
- njoy_kws: keyword argument to pass to “tape.get_pendf()”.
- **kwargskeyword argument to pass to “tape.get_ace()” plus keywords
to pass to “endf6_perturb_worker”.
- Returns:
- A dictionary of endf/pendf file or ace files depending on to_ace.
Notes
Note
ACE file temperature. Two options are implemented: - Generation of pendf file at 0K and broadening to the
required temperature when ACE file is created.
Generation of pendf file at temperature and broadening not performed when ACE is created. This approach takes into account implicit effect.
Examples
Example to produce and apply perturbations to Pu-239 xs and nubar. >>> tape = sandy.get_endf6_file(“jeff_33”, “xs”, 942390) >>> smps = tape.get_perturbations(2, njoy_kws=dict(err=1, chi=False, mubar=False, errorr33_kws=dict(mt=[2, 4, 18]),), smp_kws=dict(seed31=1, seed33=3))
Let’s apply both nubar and xs perturbations, then only nubar and then only xs. >>> outs_31_33 = tape.apply_perturbations(smps, njoy_kws=dict(err=1), processes=1) >>> outs_31 = tape.apply_perturbations({31: smps[31]}, njoy_kws=dict(err=1), processes=1) >>> outs_33 = tape.apply_perturbations({33: smps[33]}, njoy_kws=dict(err=1), processes=1)
Check that files are different for different samples. >>> for i in range(2): … assert(outs_33[i][“endf6”].data == tape.data) … assert(outs_31[i][“endf6”].data != tape.data) … assert(outs_31[i][“endf6”].data == outs_31_33[i][“endf6”].data) … assert(outs_33[i][“pendf”].data != outs_31[i][“pendf”].data) … assert(outs_33[i][“pendf”].data == outs_31_33[i][“pendf”].data)
Check that method is consistent only nubar, only xs or both nubar and xs are perturbed. >>> assert outs_33[0][“pendf”].data != outs_33[1][“pendf”].data >>> assert outs_33[0][“endf6”].data == outs_33[1][“endf6”].data >>> assert outs_31[0][“pendf”].data == outs_31[1][“pendf”].data >>> assert outs_31[0][“endf6”].data != outs_31[1][“endf6”].data
Check that redundant nubar is also perturbed. >>> nu0 = sandy.Xs.from_endf6(outs_31[0][“endf6”].filter_by(listmt=[452, 455, 456])) >>> nu1 = sandy.Xs.from_endf6(outs_31[1][“endf6”].filter_by(listmt=[452, 455, 456])) >>> assert not nu0.data[9437, 452].equals(nu1.data[9437, 452]) >>> assert nu0.data[9437, 455].equals(nu1.data[9437, 455]) >>> assert not nu0.data[9437, 456].equals(nu1.data[9437, 456])
Check that redundant and partial cross sections are correctly perturbed. >>> xs0 = sandy.Xs.from_endf6(outs_33[0][“pendf”].filter_by(listmf=[3])) >>> xs1 = sandy.Xs.from_endf6(outs_33[1][“pendf”].filter_by(listmf=[3])) >>> assert not xs0.data[9437, 1].equals(xs1.data[9437, 1]) >>> assert not xs0.data[9437, 2].equals(xs1.data[9437, 2]) >>> assert not xs0.data[9437, 4].equals(xs1.data[9437, 4]) >>> assert not xs0.data[9437, 18].equals(xs1.data[9437, 18]) >>> assert not xs0.data[9437, 51].equals(xs1.data[9437, 51]) >>> assert xs0.data[9437, 16].equals(xs1.data[9437, 16]) >>> assert xs0.data[9437, 102].equals(xs1.data[9437, 102]) >>> assert xs0.data[9437, 103].equals(xs1.data[9437, 103]) >>> assert xs0.data[9437, 107].equals(xs1.data[9437, 107])
Check that ENDF6 and PENDF output filenames are correct >>> endf6 = sandy.get_endf6_file(‘jeff_33’, ‘xs’, 10010) >>> smps = endf6.get_perturbations(2, njoy_kws=dict(err=0.1)) >>> outs = endf6.apply_perturbations(smps, to_file=True) >>> assert outs[0][“endf6”] == ‘1001_0.endf6’ and os.path.isfile(‘1001_0.endf6’) >>> assert outs[0][“pendf”] == ‘1001_0.pendf’ and os.path.isfile(‘1001_0.endf6’) >>> assert outs[1][“endf6”] == ‘1001_1.endf6’ and os.path.isfile(‘1001_1.endf6’) >>> assert outs[1][“pendf”] == ‘1001_1.pendf’ and os.path.isfile(‘1001_1.pendf’)
Check that ACE output filenames are correct >>> outs = endf6.apply_perturbations(smps, to_file=True, to_ace=True, ace_kws=dict(err=1, temperature=300, purr=False, heatr=False, thermr=False, gaspr=False)) >>> assert outs[0][“ace”] == ‘1001_0.03c’ and os.path.isfile(‘1001_0.03c’) >>> assert outs[0][“xsdir”] == ‘1001_0.03c.xsd’ and os.path.isfile(‘1001_0.03c.xsd’) >>> assert outs[1][“ace”] == ‘1001_1.03c’ and os.path.isfile(‘1001_1.03c’) >>> assert outs[1][“xsdir”] == ‘1001_1.03c.xsd’ and os.path.isfile(‘1001_1.03c.xsd’)
- get_ace(temperature=0, err=0.001, minimal_processing=False, verbose=False, **kwargs)
- Parameters:
- errTYPE, optional
reconstruction tolerance for RECONR, BROADR and THERMR. The default is 0.001.
- minimal_processing: `bool`, optional
deactivate modules THERMR, GASPR, HEATR, PURR and UNRESR. The default is False.
- temperaturefloat, optional
temperature of the cross sections in K. If not given, stop the processing after RECONR (before BROADR). The default is 0.
- verbosebool, optional
flag to print NJOY input file to screen before running the executable. The default is False.
- get_errorr(temperature=0, err=0.001, minimal_processing=False, verbose=False, **kwargs)
- Parameters:
- errTYPE, optional
reconstruction tolerance for RECONR, BROADR and THERMR. The default is 0.001.
- minimal_processing: `bool`, optional
deactivate modules THERMR, GASPR, HEATR, PURR and UNRESR. The default is False.
- temperaturefloat, optional
temperature of the cross sections in K. If not given, stop the processing after RECONR (before BROADR). The default is 0.
- verbosebool, optional
flag to print NJOY input file to screen before running the executable. The default is False.
- get_gendf(temperature=0, err=0.001, minimal_processing=False, verbose=False, **kwargs)
- Parameters:
- errTYPE, optional
reconstruction tolerance for RECONR, BROADR and THERMR. The default is 0.001.
- minimal_processing: `bool`, optional
deactivate modules THERMR, GASPR, HEATR, PURR and UNRESR. The default is False.
- temperaturefloat, optional
temperature of the cross sections in K. If not given, stop the processing after RECONR (before BROADR). The default is 0.
- verbosebool, optional
flag to print NJOY input file to screen before running the executable. The default is False.
- get_id(method='nndc')
Extract ID for a given MAT for a ENDF-6 file.
- Parameters:
- methodstr, optional
Methods adopted to produce the ID. The default is “nndc”. - if method=’aleph’ the ID is the ZAM identifier - else, the ID is the ZA identifier according to the NNDC rules
- Returns:
- IDint
ID of the ENDF-6 file.
Notes
Note
a warning is raised if more than one MAT is found. Only the ID corresponding to the lowest MAT will be returned.
Examples
Extract ID for H1 file using NNDC and ALEPH methods >>> tape = sandy.get_endf6_file(“jeff_33”, “xs”, 10010) >>> assert tape.get_id() == 1001 >>> assert tape.get_id(method=”aleph”) == 10010
Extract ID for Am242m file using NNDC and ALEPH methods >>> tape2 = sandy.get_endf6_file(“jeff_33”, “xs”, 952421) >>> assert tape2.get_id() == 95642 >>> assert tape2.get_id(method=”ALEPH”) == 952421
>>> assert tape.merge(tape2).get_id() == 1001 >>> assert tape2.merge(tape).get_id() == 1001
- get_nsub()
Determine ENDF-6 sub-library type by reading flag “NSUB” of first MAT in file.
- Returns:
- int
NSUB value
Examples
assert sandy.get_endf6_file(“jeff_33”, “xs”, 10010).get_nsub() == “neutron” assert sandy.get_endf6_file(“jeff_33”, “xs”, 10010).get_nsub().get_pendf(err=1).get_nsub() == “neutron” assert sandy.get_endf6_file(“jeff_33”, “nfpy”, 942410).get_nsub() == “nfpy” assert sandy.get_endf6_file(“jeff_33”, “decay”, 942410).get_nsub() == “decay” assert sandy.get_endf6_file(“jeff_33”, “dxs”, 26000).get_nsub() == “neutron” assert sandy.get_endf6_file(“proton”, “dxs”, 26000).get_nsub() == “proton”
- get_pendf(temperature=0, err=0.001, minimal_processing=False, verbose=False, **kwargs)
- Parameters:
- errTYPE, optional
reconstruction tolerance for RECONR, BROADR and THERMR. The default is 0.001.
- minimal_processing: `bool`, optional
deactivate modules THERMR, GASPR, HEATR, PURR and UNRESR. The default is False.
- temperaturefloat, optional
temperature of the cross sections in K. If not given, stop the processing after RECONR (before BROADR). The default is 0.
- verbosebool, optional
flag to print NJOY input file to screen before running the executable. The default is False.
- get_perturbations(nsmp, to_excel=None, njoy_kws={}, smp_kws={}, **kwargs)
Construct multivariate distributions with a unit vector for mean and with relative covariances taken from the evaluated files processed with the NJOY module ERRORR.
Perturbation factors are sampled with the same multigroup structure of the covariance matrix and are returned by nuclear datatype as a dict of pd.Dataframe instances .
- Parameters:
- nsmpTYPE
DESCRIPTION.
- to_excelTYPE, optional
DESCRIPTION. The default is None.
- njoy_kwsTYPE, optional
DESCRIPTION. The default is {}.
- smp_kwsTYPE, optional
DESCRIPTION. The default is {}.
- **kwargsTYPE
DESCRIPTION.
- Returns:
- smpTYPE
DESCRIPTION.
Examples
Generate a couple of samples from the H1 file of JEFF-3.3. >>> njoy_kws = dict(err=1, errorr_kws=dict(mt=102)) >>> tape = sandy.get_endf6_file(“jeff_33”, “xs”, 10010) >>> smps = tape.get_perturbations(nsmp=2, njoy_kws=njoy_kws) >>> assert len(smps) == 1 >>> assert isinstance(smps[33], sandy.Samples) >>> assert (smps[33].data.index.get_level_values(“MT”) == 102).all()
- get_records()
Extract MAT, MF and MT combinations avaialbel in the file and report it in tabulated format.
- Returns:
- dfpd.DataFrame
Dataframe with MAT, MF and MT as columns.
Examples
Short test for hydrogen >>> sandy.get_endf6_file(“jeff_33”, “xs”, 10010).get_records()
MAT MF MT
0 125 1 451 1 125 2 151 2 125 3 1 3 125 3 2 4 125 3 102 5 125 4 2 6 125 6 102 7 125 33 1 8 125 33 2 9 125 33 102
- read_section(mat, mf, mt, raise_error=True)
Parse MAT/MF/MT section.
- Parameters:
- `mat`int
MAT number
- `mf`int
MF number
- `mt`int
MT number
- Returns:
- dict
- update_intro(**kwargs)
Method to update MF1/MT451 of each MAT based on the file content (concistency is enforced) and user-given keyword arguments.
- Parameters:
- **kwargsdict
dictionary of elements to be modified in section MF1/MT451 (it applies to all MAT numbers).
- Returns:
Endf6()
Endf6()
with updated MF1/MT451.
Examples
Check how many lines of description and how many sections are recorded in a file. >>> tape = sandy.get_endf6_file(“jeff_33”, “xs”, 10010) >>> intro = tape.read_section(125, 1, 451) >>> assert len(intro[“DESCRIPTION”]) == 87 >>> assert len(intro[“SECTIONS”]) == 10
By removing sections in the Endf6 instance, the recorded number of sections does not change. >>> tape2 = tape.delete_section(125, 33, 1).delete_section(125, 33, 2) >>> intro = tape2.read_section(125, 1, 451) >>> assert len(intro[“DESCRIPTION”]) == 87 >>> assert len(intro[“SECTIONS”]) == 10
Running updated intro updates the recorded number of sections. >>> tape2 = tape.delete_section(125, 33, 1).delete_section(125, 33, 2).update_intro() >>> intro = tape2.read_section(125, 1, 451) >>> assert len(intro[“DESCRIPTION”]) == 87 >>> assert len(intro[“SECTIONS”]) == 8
It can also be used to update the lines of description. >>> intro = tape2.update_intro(**dict(DESCRIPTION=[” new description”])).read_section(125, 1, 451) >>> print(sandy.write_mf1(intro))
1001.00000 9.991673-1 0 0 2 5 125 1451 1 0.00000000 0.00000000 0 0 0 6 125 1451 2 1.00000000 20000000.0 3 0 10 3 125 1451 3 0.00000000 0.00000000 0 0 1 8 125 1451 4 new description 125 1451 5
1 451 13 0 125 1451 6 2 151 4 0 125 1451 7 3 1 35 0 125 1451 8 3 2 35 0 125 1451 9 3 102 35 0 125 1451 10 4 2 196 0 125 1451 11 6 102 201 0 125 1451 12
33 102 21 0 125 1451 13
- sandy.core.endf6.get_endf6_file(library, kind, zam, to_file=False)
Given a library and a nuclide import the corresponding ENDF-6 nuclear data file directly from internet.
- Parameters:
- librarystr
nuclear data library. Available libraries are: for ‘xs’:
‘endfb_71’
‘endfb_80’
‘irdff_2’
‘jeff_32’
‘jeff_33’
‘jendl_40u’
- for ‘nfpy’:
‘endfb_71’
‘jeff_311’
‘jeff_33’
‘endfb_80’
‘jendl_40u’
- for ‘decay’:
‘endfb_71’
‘jeff_311’
‘jeff_33’
‘endfb_80’
- for ‘tsl’: (read the note)
‘endfb_71’
‘jeff_33’
‘endfb_80’
‘jendl_40u’
- for ‘dxs’:
‘jeff_33’
‘proton’
- kindstr
- nuclear data type:
‘xs’ is a standard neutron-induced nuclear data file
‘nfpy’ is a Neutron-Induced Fission Product Yields nuclear data file
‘decay’ is a Radioactive Decay Data nuclear data file
‘dxs’ is displacement cross section data file
‘tsl’ is a Thermal Neutron Scattering Data file
- zamint or ‘all’ or iterable
- zam = ‘int’ (individual nuclides) or iterable (group of nuclides)
- ZAM nuclide identifier $Z times 10000 + A times 10 + M$ where:
$Z$ is the charge number
$A$ is the mass number
$M$ is the metastate level (0=ground, 1=1st level)
- zam = ‘all’
We obtain the information of all the library. This option is not available for ‘xs’
- Returns:
- Endf6
Endf6 object with ENDF-6 data for specified library and nuclide.
- Raises:
- ValueError
if library is not among available selection.
- ValueError
if when you select ‘xs’, you select zam = ‘all’
Notes
Note
In the kind=’tls’ option, instead of the zam, integers are used. If you need help, the get_tsl_index function contains all the necessary information for the correct choice of these integers.
Examples
Import hydrogen file from JEFF-3.3. >>> tape = sandy.get_endf6_file(“jeff_33”, ‘xs’, 10010) >>> assert type(tape) is sandy.Endf6
Import hydrogen file from ENDF/B-VII.1. >>> tape = sandy.get_endf6_file(“endfb_71”, ‘xs’, 10010) >>> assert type(tape) is sandy.Endf6
Import hydrogen file from ENDF/B-VIII.0. >>> tape = sandy.get_endf6_file(“endfb_80”, ‘xs’, 10010) >>> assert type(tape) is sandy.Endf6
Import hydrogen file from JENDL-4.0u >>> tape = sandy.get_endf6_file(“jendl_40u”, ‘xs’, 10010) >>> assert type(tape) is sandy.Endf6
Import Neutron-Induced Fission Product Yields for Th-227 from ENDF/B-VII.1. >>> tape = sandy.get_endf6_file(“endfb_71”, ‘nfpy’, 902270) >>> assert type(tape) is sandy.Endf6
Import Neutron-Induced Fission Product Yields for Th-227 from ENDF/B-VIII.0 >>> tape = sandy.get_endf6_file(“endfb_80”, ‘nfpy’, 902270) >>> assert type(tape) is sandy.Endf6
Import Neutron-Induced Fission Product Yields for Th-227 from JENDL-4.0u >>> tape = sandy.get_endf6_file(“jendl_40u”, ‘nfpy’, 902270) >>> assert type(tape) is sandy.Endf6
Import Neutron-Induced Fission Product Yields for Th-232 from JEFF-3.1.1 >>> tape = sandy.get_endf6_file(“jeff_311”, ‘nfpy’, 902320) >>> assert type(tape) is sandy.Endf6
Import Neutron-Induced Fission Product Yields for Th-232 from JEFF-3.3 >>> tape = sandy.get_endf6_file(“jeff_33”, ‘nfpy’, 902320) >>> assert type(tape) is sandy.Endf6
Import Radioactive Decay Data for H-1 from JEFF-3.1.1 >>> tape = sandy.get_endf6_file(“jeff_311”, ‘decay’, 10010) >>> assert type(tape) is sandy.Endf6
Import Radioactive Decay Data for H-1 from JEFF-3.3 >>> tape = sandy.get_endf6_file(“jeff_33”, ‘decay’, 10010) >>> assert type(tape) is sandy.Endf6
Import Radioactive Decay Data for H-1 from ENDF/B-VII.1. >>> tape = sandy.get_endf6_file(“endfb_71”, ‘decay’, 10010) >>> assert type(tape) is sandy.Endf6
Import Radioactive Decay Data for H-1 from ENDF/B-VIII.0. >>> tape = sandy.get_endf6_file(“endfb_80”, ‘decay’, 10010) >>> assert type(tape) is sandy.Endf6
Import all Neutron-Induced Fission Product Yields from ENDF/B-VII.1. >>> tape = sandy.get_endf6_file(“endfb_71”, ‘nfpy’, ‘all’) >>> assert type(tape) is sandy.Endf6
Import a list of Decay Data for JEFF-3.3. >>> tape = sandy.get_endf6_file(“jeff_33”, ‘decay’, [380900, 551370, 541350]) >>> assert type(tape) is sandy.Endf6
Thermal Neutron Scattering Data from ENDF/B-VII.1. >>> tape = sandy.get_endf6_file(“endfb_71”, ‘tsl’, [1, 2, 3]) >>> assert type(tape) is sandy.Endf6
Thermal Neutron Scattering Data from ENDF/B-VIII.0. >>> tape = sandy.get_endf6_file(“endfb_80”, ‘tsl’, [1, 2, 3]) >>> assert type(tape) is sandy.Endf6
Thermal Neutron Scattering Data from JEFF-3.3. >>> tape = sandy.get_endf6_file(“jeff_33”, ‘tsl’, [1, 2, 3]) >>> assert type(tape) is sandy.Endf6
Thermal Neutron Scattering Data from JENDL-4.0u >>> tape = sandy.get_endf6_file(“jendl_40u”, ‘tsl’, [1, 2, 3]) >>> assert type(tape) is sandy.Endf6
Import natural Fe for IRDFF-II >>> tape = sandy.get_endf6_file(“irdff_2”, “xs”, 260000) >>> assert type(tape) is sandy.Endf6
- sandy.core.endf6.get_tsl_index(library)
Obtain the index information available in the library web page.
- Parameters:
- librarystr
nuclear data library. Available libraries are: for ‘tsl’
‘endfb_71’
‘jeff_33’
‘endfb_80’
‘jendl_40u
‘irdff_ii
- Raises:
- ValueError
if library is not among available selection.
- Example
>>> sandy.endf6.get_tsl_index("jendl_40u") Lib: JENDL-4.0 Library: JENDL-4.0 Japanese evaluated nuclear data library, 2010 Sub-library: NSUB=12 Thermal Neutron Scattering Data
KEY Material Lab. Date Authors
1 1-H(H2O) LANL EVAL-apr93 MACFARLANE 20.MeV tsl_0001_h(h2o).zip 412Kb
2 1-Para-H LANL EVAL-APR93 MacFarlane 20.MeV tsl_0002_para-H.zip 91Kb
3 1-Ortho-H LANL EVAL-APR93 MacFarlane 20.MeV tsl_0003_ortho-H.zip 96Kb
7 1-H(ZrH) LANL EVAL-apr93 MACFARLANE 20.MeV tsl_0007_h(zrh).zip 448Kb
11 1-D(D2O) GA EVAL-DEC69 KOPPEL,HOUSTON 20.MeV tsl_0011_D(D2O).zip 235Kb
12 1-Para-D LANL EVAL-APR93 MacFarlane 20.MeV tsl_0012_para-d.zip 92Kb
13 1-Ortho-D LANL EVAL-APR93 MacFarlane 20.MeV tsl_0013_ortho-d.zip 93Kb
26 4-Be-metal LANL EVAL-apr93 MACFARLANE 20.MeV tsl_0026_bemetal.zip 419Kb
27 4-BeO LANL EVAL-apr93 MACFARLANE 20.MeV tsl_0027_beo.zip 483Kb
31 6-Graphite LANL EVAL-apr93 MACFARLANE 20.MeV tsl_0031_graphite.zip 397Kb
33 6-l-CH4 LANL EVAL-APR93 MacFarlane 20.MeV tsl_0033_l-ch4.zip 50Kb
34 6-s-CH4 LANL EVAL-APR93 MacFarlane 20.MeV tsl_0034_s-ch4.zip 42Kb
37 6-H(CH2) GA EVAL-DEC69 KOPPEL,HOUSTON,SPREVAK 20.MeV tsl_0037_H(CH2).zip 72Kb
40 6-BENZINE GA EVAL-DEC69 KOPPEL,HOUSTON,BORGONOVI 20.MeV tsl_0040_BENZINE.zip 236Kb
58 40-Zr(ZrH) LANL EVAL-apr93 MACFARLANE 20.MeV tsl_0058_zr(zrh).zip 201Kb
- Total: Materials:15 Size:11Mb Compressed:4Mb