fudge.processing.resonances package

Submodules

fudge.processing.resonances.getCoulombWavefunctions module

fudge.processing.resonances.getCoulombWavefunctions.coulombNormalizationFactor(L, eta)[source]

Coulomb wavefunction normalization factor (see DLMF Eq. 33.2.5 or Abramowitz & Stegun Eq. 14.1.7), C_\ell(\eta)={2^\ell e^{2\pi\eta}|\Gamma(\ell+1+i\eta)|}/{(2\ell+1)!}

Because we’re working with numpy arrays of \eta (for parallelization ala’ Caleb’s tricks), the logic of numpy masked arrays is at work here to handle different regimes of eta.

Parameters:
  • L (int) – the L to evaluate at
  • eta (numpy.array(type=float)) – array of \eta values
Returns:

array of normalizations

Return type:

numpy.array(type=float)

fudge.processing.resonances.getCoulombWavefunctions.coulombPenetrationFactor(L, rho, eta)[source]

Compute the Coulomb penetrability, P_\ell, P_\ell(\rho,\eta)={\rho}/{(A_\ell(\rho,\eta))^2} where (A_\ell(\rho,\eta))^2=(G_\ell(\rho,\eta))^2+(F_\ell(\rho,\eta))^2

Because we’re working with numpy arrays of \eta (for parallelization ala’ Caleb’s tricks), the logic of numpy masked arrays is at work here to handle different regimes of eta.

Here we use an external subroutine coulfg2 from Thompson et al, converted to c and wrapped

Parameters:
  • L (int) – the L to evaluate at
  • rho (numpy.array(type=float)) – array of \rho values
  • eta (numpy.array(type=float)) – array of \eta values
Returns:

array of penetrabilities

Return type:

numpy.array(type=float)

fudge.processing.resonances.getCoulombWavefunctions.coulombPhi(L, rho, eta)[source]

Compute the Coulomb phase, \varphi_\ell, \varphi_\ell(\rho,\eta)=\arg(G_\ell(\rho,\eta)+iF_\ell(\rho,\eta))

Because we’re working with numpy arrays of \eta (for parallelization ala’ Caleb’s tricks), the logic of numpy masked arrays is at work here to handle different regimes of eta.

Here we use an external subroutine ‘coulfg2’ from Thompson et al, converted to c and wrapped

Parameters:
  • L (int) – the L to evaluate at
  • rho (numpy.array(type=float)) – array of \rho values
  • eta (numpy.array(type=float)) – array of \eta values
Returns:

numpy.array(type=float) of phases

fudge.processing.resonances.getCoulombWavefunctions.coulombShiftFactor(L, rho, eta)[source]

Compute the Coulomb shift, S_\ell, S_\ell(\rho,\eta)=({\rho}/{A_\ell(\rho,\eta)})({\partial A_\ell(\rho,\eta)}/{\partial\rho}) where (A_\ell(\rho,\eta))^2=(G_\ell(\rho,\eta))^2+(F_\ell(\rho,\eta))^2

Because we’re working with numpy arrays of \eta (for parallelization ala’ Caleb’s tricks), the logic of numpy masked arrays is at work here to handle different regimes of eta.

Here we use an external subroutine ‘coulfg2’ from Thompson et al, converted to c and wrapped

Parameters:
  • L (int) – the L to evaluate at
  • rho (numpy.array(type=float)) – array of \rho values
  • eta (numpy.array(type=float)) – array of \eta values
Returns:

array of shifts

Return type:

numpy.array(type=float)

fudge.processing.resonances.getCoulombWavefunctions.digamma(n)[source]

Simple digamma function, tested against A&S tables in chapter 6, up to n=20 :param n: an integer :return:

fudge.processing.resonances.getCoulombWavefunctions.getCoulombWavefunctions(rho, eta, L)[source]
fudge.processing.resonances.getCoulombWavefunctions.test_getCoulombWavefunctions()[source]

fudge.processing.resonances.getScatteringMatrices module

fudge.processing.resonances.getScatteringMatrices.getScatteringMatrices(E, Eres, captureWidth, widths, penetrabilities)[source]

fudge.processing.resonances.reconstructResonances module

class fudge.processing.resonances.reconstructResonances.ChannelDesignator(l, J, reaction, index, s, gfact=None, particleA=None, particleB=None, Xi=0.0, isElastic=None, channelClass=None, useRelativistic=False, eliminated=False)[source]

Bases: object

J
Xi
channelClass
eliminated
gfact
has_same_J(other)[source]
index
isElastic
is_open(Ein)[source]
l
particleA
particleB
reaction
s
useRelativistic
class fudge.processing.resonances.reconstructResonances.MLBWcrossSection(reactionSuite, sectionIndex=None, enableAngDists=False, verbose=True, **kw)[source]

Bases: fudge.processing.resonances.reconstructResonances.RRBaseClass

Given a resonance region in MLBW format, create a class with all data required to reconstruct

Only the elastic channel differs from SLBW

getChannelConstantsBc()[source]

For ENDF’s MLBW, should be B_c = S_\ell(|E_\lambda|) where \ell is the channel angular momentum and \lambda is the resonances index for the channel

getCrossSection(E, **kwargs)[source]
getScatteringMatrixU(Ein, useTabulatedScatteringRadius=True, enableExtraCoulombPhase=False)[source]

Compute the scattering matrix. We could have used the generic U function in the base class, but Froehner has “simplifications” that we took advantage of here (that and I don’t know what the R matrix is exactly in the case of MLBW).

class fudge.processing.resonances.reconstructResonances.RMatrixLimitedcrossSection(reactionSuite, sectionIndex=None, enableAngDists=False, verbose=True, **kw)[source]

Bases: fudge.processing.resonances.reconstructResonances.RRBaseClass

extended Reich_Moore (LRF=7 in ENDF) Here, resonances are sorted primarily by J: within each ‘spin group’, total J is conserved One or more competitive channels may be used in this case. Also, each resonance may have contributions from multiple l-waves

eta(Ex, pA, pB)[source]

eta, the Sommerfeld parameter, given by \eta = Z_A Z_B m_{red} \alpha / ( \hbar c k )

for competitive channels with 2 charged particles, parameter eta is used to find penetrability. eta is given in eq D.79 of ENDF manual and $e^2$ is the fine structure constant $lpha$ and $m_{red}$ is the reduced mass. eta is dimensionless.

Parameters:
  • Ex – The incident energy
  • pA – particle A
  • pB – particle B
Returns:

the Sommerfeld parameter [dimensionless]

getAPByChannel(c, trueOrEffective='true')[source]
getChannelConstantsBc()[source]

For ENDF’s Reich-Moore, should be B_c = -\ell where \ell is the channel angular momentum, but the ENDF manual says nothing about it.

There is a per-channel parameter BCH that we will interpret as B_c

getChannelScatteringRadiiEffective()[source]
getChannelScatteringRadiiTrue()[source]
getCrossSection(E, **kwargs)[source]
getEiPhis(Ein, useTabulatedScatteringRadius=None, enableExtraCoulombPhase=False)[source]

The phase factor for the collision matrix, \Omega_c=e^{\omega_c-\varphi_c}

Parameters:
  • Ein (numpy.array(type=float)) –
  • useTabulatedScatteringRadius (bool) –
  • enableExtraCoulombPhase (bool) –
Returns:

Return type:

getL0Matrix(Ein)[source]

Get the L^0 matrix of Froehner, {\bf L^0}_{cc'} = \delta_{cc'} (L_c-B_c) where L_c = S_c + i P_c

getLMax(maxLmax=10)[source]

LMax is determined by the behavior of the Blatt-Biedenharn Zbar coefficients. Inside each one, there is a Racah coefficient and a Clebsh-Gordon coefficient. The CG coefficient looks like this:

( l1 l2 L ) ( 0 0 0 )

So, this means two things. First, the CG coeff (and hence Zbar) will be zero if l1+l2+L=odd. Second, the maximum value of L will be l1max+l2max. Hence, Lmax=2*lmax.

getRMatrix(Ein)[source]

The R matrix in the Reich-Moore approximation is

System Message: WARNING/2 (R_{cc'}=\sum_\lambda{\gamma_{\lambda c}\gamma_{\lambda c'}/({E_\lambda - E - i\Gamma_{\lambda\gamma}/2}))

latex exited with error [stdout] This is pdfTeX, Version 3.14159265-2.6-1.40.16 (TeX Live 2015) (preloaded format=latex) restricted \write18 enabled. entering extended mode (./math.tex LaTeX2e <2015/01/01> Babel <3.9l> and hyphenation patterns for 79 languages loaded. (/usr/local/texlive/2015/texmf-dist/tex/latex/base/article.cls Document Class: article 2014/09/29 v1.4h Standard LaTeX document class (/usr/local/texlive/2015/texmf-dist/tex/latex/base/size12.clo)) (/usr/local/texlive/2015/texmf-dist/tex/latex/base/inputenc.sty (/usr/local/texlive/2015/texmf-dist/tex/latex/ucs/utf8x.def)) (/usr/local/texlive/2015/texmf-dist/tex/latex/ucs/ucs.sty (/usr/local/texlive/2015/texmf-dist/tex/latex/ucs/data/uni-global.def)) (/usr/local/texlive/2015/texmf-dist/tex/latex/amsmath/amsmath.sty For additional information on amsmath, use the `?' option. (/usr/local/texlive/2015/texmf-dist/tex/latex/amsmath/amstext.sty (/usr/local/texlive/2015/texmf-dist/tex/latex/amsmath/amsgen.sty)) (/usr/local/texlive/2015/texmf-dist/tex/latex/amsmath/amsbsy.sty) (/usr/local/texlive/2015/texmf-dist/tex/latex/amsmath/amsopn.sty)) (/usr/local/texlive/2015/texmf-dist/tex/latex/amscls/amsthm.sty) (/usr/local/texlive/2015/texmf-dist/tex/latex/amsfonts/amssymb.sty (/usr/local/texlive/2015/texmf-dist/tex/latex/amsfonts/amsfonts.sty)) (/usr/local/texlive/2015/texmf-dist/tex/latex/anyfontsize/anyfontsize.sty) (/usr/local/texlive/2015/texmf-dist/tex/latex/tools/bm.sty) No file math.aux. (/usr/local/texlive/2015/texmf-dist/tex/latex/ucs/ucsencs.def) (/usr/local/texlive/2015/texmf-dist/tex/latex/amsfonts/umsa.fd) (/usr/local/texlive/2015/texmf-dist/tex/latex/amsfonts/umsb.fd) ! Missing } inserted. <inserted text> } l.13 ...\lambda - E - i\Gamma_{\lambda\gamma}/2})$ [1] (./math.aux) ) (see the transcript file for additional information) Output written on math.dvi (1 page, 572 bytes). Transcript written on math.log.

Parameters:Ein (numpy.array(type=float)) –
Returns:
Return type:
getScatteringMatrixT(Ein, useTabulatedScatteringRadius=True, enableExtraCoulombPhase=False)[source]
Parameters:
  • Ein (numpy.array(type=float)) –
  • useTabulatedScatteringRadius (bool) –
  • enableExtraCoulombPhase (bool) –
Returns:

Return type:

isElastic(reactionDesignator)[source]
k_competitive(Ex, pA, pB)[source]

Calculate k for any 2-body output channel.

Note that if pA and pB are target and neutron, this reduces to self.k(E) as defined above in the resonanceReconstructionBaseClass

Parameters:
  • Ex – incident energy - Xi (Xi is the lab frame reaction threshold) in eV
  • pA – particle A
  • pB – particle B
Returns:

k in b**-1/2

omega(eta, L)[source]
penetrationFactorByChannel(c, Ein)[source]
phiByChannel(c, Ein)[source]
resetResonanceParametersByChannel(multipleSScheme='ENDF', useReichMooreApproximation=False, Ein=None)[source]
rho(Ein, c)[source]

Compute \rho_c(E) = a_c * k_c(E), using the true scattering radius. ENDF uses it for calculating shift and penetrabilities.

Parameters:
  • Ein – incident energy in the lab frame (shifted by a threshold, if appropriate)
  • c – the channel designator
Returns:

the value of rho (dimensionless)

rhohat(Ein, c)[source]

Compute \hat{\rho}_c(E) = a_c * k_c(E), using the effective scattering radius ENDF uses it for calculating the phase (but in truth, there should be no effective scattering radius). (Caleb uses self.k below, but I think it should be self.k_competitive for the sake of consistency)

Parameters:
  • Ein – incident energy in the lab frame (shifted by a threshold, if appropriate)
  • c – the channel designator
Returns:

the value of rho (dimensionless)

setResonanceParametersByChannel(multipleSScheme='ENDF', useReichMooreApproximation=False, Ein=None, warnOnly=False)[source]

Reorganize member data into channels :param multipleSScheme: ignored, kept so has same signature as overridden function :param useReichMooreApproximation: ignored, kept so has same signature as overridden function :return:

shiftFactorByChannel(c, Ein)[source]
class fudge.processing.resonances.reconstructResonances.RMcrossSection(reactionSuite, sectionIndex=None, enableAngDists=False, verbose=False, **kw)[source]

Bases: fudge.processing.resonances.reconstructResonances.RRBaseClass

simplified Reich_Moore (LRF=3 in ENDF) More complex than Breit-Wigner approximations, but we can still use same __init__

getChannelConstantsBc()[source]

For ENDF’s Reich-Moore, should be

..math::
B_c = -ell

where $ell$ is the channel angular momentum

what is not clearly stated in the ENDF manual is that Bc = 0 for this case

getCrossSection(E, **kwargs)[source]
getL0Matrix(Ein, Bc=None)[source]

Get the L0 matrix of Froehner:

..math:

{f L^0}_{cc'} = \delta_{cc'} (L_c-B_c)

where

..math:

L_c = S_c + i P_c

But…. ENDF’s RM formulation uses a shift of 0 and a Bc of 0!!!

getRMatrix(Ein)[source]
The R matrix in the Reich-Moore approximation, :math:`R_{cc’} = sum_lambda

rac{gamma_{lambda c}gamma_{lambda c’}{E_lambda - E - iGamma_{lambdagamma}/2}`

class fudge.processing.resonances.reconstructResonances.RRBaseClass(reactionSuite, sectionIndex=None, lowerBound=None, upperBound=None, RR=None, **kw)[source]

Bases: fudge.processing.resonances.reconstructResonances.resonanceReconstructionBaseClass

generateEnergyGrid()[source]

Create an initial energy grid by merging a rough mesh for the entire region (~10 points / decade) with a denser grid around each resonance. For the denser grid, multiply the total resonance width by the ‘resonancePos’ array defined below.

getAMatrix(Ein)[source]
getAPByChannel(c, trueOrEffective='true')[source]

get the channel radius, rho. If L is specified try to get L-dependent value

getAngularDistribution(E, **kwargs)[source]
getAverageQuantities(nBins=1, binScheme='linspacing', computeUncertainty=False)[source]

Computes average widths and level spacings from the set of resonance parameters in self, on a per-channel basis

The averages are computed in equal lethargy bins starting at the lowest resonance energy in a sequence up to the upperBound of the resonance region. I tried to keep on average 10 resonances/logrithmic bin so I can get a reasonable average.

Parameters:
  • nBins – the number of resonances per logrithmic bin to aim for, on average
  • binScheme – the scheme to use to determine the binning
  • computeUncertainty – toggle the calculation of the uncertainty of the quantities
Returns:

a dictionary of results, sorted by channel

getBackgroundRMatrix(Ein, Emin, Emax, gamWidth, pole_strength, Rinf=None)[source]

Froehner’s background R matrix. This could be vectorized for speed. Not sure if it would be useful though.

getChannelConstantsBc()[source]
getEiPhis(Ein, useTabulatedScatteringRadius=True, enableExtraCoulombPhase=False)[source]

The phase factor for the collision matrix:

..math::
Omega_c = e^{-

arphi_c}

getKMatrix(Ein)[source]
getL0Matrix(Ein)[source]

Get the L0 matrix of Froehner:

..math::
{f L^0}_{cc’} = delta_{cc’} (L_c-B_c)

where

..math::
L_c = S_c + i P_c
getLMax(maxLmax=10)[source]

LMax is determined by the behavior of the Blatt-Biedenharn Zbar coefficients. Inside each one, there is a Racah coefficient and a Clebsh-Gordon coefficient. The CG coefficient looks like this:

( l1 l2 L ) ( 0 0 0 )

So, this means two things. First, the CG coeff (and hence Zbar) will be zero if l1+l2+L=odd. Second, the maximum value of L will be l1max+l2max. Hence, Lmax=2*lmax.

getParticleParities(rxn)[source]

Wrapper around getParticleSpinParities to extract just the parities

Parameters:rxn – the reaction string for this resonanceReaction (or equivalent). We’ll process this string to out what the light particle and heavy “residual” is and then look up their JPi.
Returns:a tuple: (JPi_light, JPi_heavy). JPi itself is a tuple (J,Pi). J is float (either integer or 1/2 integer) and Pi is -1 or +1.
getParticleSpinParities(rxn)[source]
6 cases:

rxn==’capture’ rxn==’elastic’ (equivalent to ‘projectile + target’) rxn==’fission’ or ‘fissionA’ or ‘fissionB’ rxn==’competative’ rxn==’something + something’ where one of somethings is photon (same as ‘capture’) rxn==’something + something’ where neither something is a photon

A word about capture channels: Because of the use of the Reich-Moore approximation, photon channels may be quasi-channels consisting of many photon channels lumped together. This makes determination of the residual nucleus spin & parity difficult since it really corresponds to a bunch of states. Here, we just compute the lowest spin consistent with the target & photon and leave it up to the code to override these values as needed.

A word about fission and competative channels: Here we don’t really know or need the spins & parities of the particles. Both are many channels lumped together, in a way that makes capture look simple. We’ll just return the elastic parameters under the assumption that they’re not needed.

Parameters:rxn – the reaction string for this resonanceReaction (or equivalent). We’ll process this string to out what the light particle and heavy “residual” is and then look up their JPi.
Returns:a tuple: (JPi_a, JPi_b). JPi itself is a tuple (J,Pi). J is float (either integer or 1/2 integer) and Pi is -1 or +1. It doesn’t matter whether particle a or b is heaviest
getParticleSpins(rxn)[source]

Wrapper around getParticleSpinParities to extract just the spins

Parameters:rxn – the reaction string for this resonanceReaction (or equivalent). We’ll process this string to out what the light particle and heavy “residual” is and then look up their JPi.
Returns:a tuple: (JPi_light, JPi_heavy). JPi itself is a tuple (J,Pi). J is float (either integer or 1/2 integer) and Pi is -1 or +1.
getPoleStrength(computeUncertainty=False)[source]

Computes the neutron pole strength from the transmission coefficients

..math::
s_c(E) =

ho_c(E)overline{Gamma_c}/2P_c = T_c(E)/4pi P_c

param computeUncertainty:
 
return:
getPorterThomasFitToWidths(Emin=0.0, Emax=None, verbose=False)[source]

Perform a channel-by-channel fit of the histogram of widths to a Porter-Thomas distribution.

Parameters:verbose
Returns:
getRMatrix(Ein)[source]
getReducedWidths(channel, Emin=None, Emax=None)[source]

Convert the list of widths to reduced widths:

redWidth = width/2*penetrability
Parameters:
  • widthList
  • channel
Returns:

getScatteringLength()[source]

Compute R’ in b**1/2, should be close to AP

The potential scattering cross section sigPot = 4 Pi (R’)^2, so we compute the potential scattering cross section at E=1e-5 eV :return:

getScatteringMatrixT(Ein, useTabulatedScatteringRadius=True, enableExtraCoulombPhase=False)[source]
getScatteringMatrixU(Ein, useTabulatedScatteringRadius=True, enableExtraCoulombPhase=False)[source]

Compute the scattering matrix using

getStrengthFunction(computeUncertainty=False)[source]
getTransmissionCoefficientsFromSumRule(resonancesPerBin=10, computeUncertainty=False)[source]
Use Moldauer’s sum rule to extract the transmission coefficients directly from the RRR tables
P.A. Moldauer Phys. Rev. Lett. 19, 1047-1048 (1967)
Parameters:computeUncertainty – toggle the calculation of the uncertainty of the transmission coefficients
Returns:a dictionary of results, sorted by channel
getWMatrix(Ein)[source]

Get the W matrix for use in computing U. W is:

..math::
{f W} = {f I} + 2i{f X}
getXMatrix(Ein)[source]

Get the X matrix for use in computing W. X is:

..math::
{f X}_{cc’} = P^{-1/2}_c ( ( {f I} - {f R}{f L^0} )^{-1}{f R} )_{cc’} P_{c’}^{-1/2}delta_{JJ’}
penetrationFactorByChannel(c, Ein)[source]
phiByChannel(c, Ein)[source]
resetResonanceParametersByChannel(multipleSScheme='ENDF', useReichMooreApproximation=False, Ein=None)[source]
rho(E, L=None)[source]

get the channel radius, rho. If L is specified try to get L-dependent value

setResonanceParametersByChannel(multipleSScheme='ENDF', useReichMooreApproximation=False, Ein=None, warnOnly=False)[source]

Reorganize member data into channels (relies heavily on groundwork in sortLandJ).

shiftFactorByChannel(c, Ein)[source]
sortLandJ()[source]

SLBW, MLBW and Reich_Moore formalisms have similar structure it’s convenient to sort their resonances by L and J This method should NOT be used for R-Matrix Limited

fudge.processing.resonances.reconstructResonances.RoundToSigFigs(x, sigfigs)[source]

Rounds the value(s) in x to the number of significant figures in sigfigs.

Restrictions: sigfigs must be an integer type and store a positive value. x must be a real value or an array like object containing only real values.

class fudge.processing.resonances.reconstructResonances.SLBWcrossSection(reactionSuite, sectionIndex=None, enableAngDists=False, verbose=True, **kw)[source]

Bases: fudge.processing.resonances.reconstructResonances.RRBaseClass

given a resonance region in SLBW format, create a class with all data required to reconstruct

Note, the resonances in the SLBW format each correspond to different “levels” and so do not interfere. This is unlike all of the other RRR formats. So, while one resonance energy and one width go with one channel and all the channels associated with one resonance energy go together to make one reaction, the reactions for each resonance are added incoherently (no interference effects at all).

getAngularDistribution(E, **kwargs)[source]
getCrossSection(E, **kwargs)[source]
getScatteringLength()[source]

Compute R’ in b**1/2, should be close to AP

The potential scattering cross section sigPot = 4 Pi (R’)^2, so we compute the potential scattering cross section at E=1e-5 eV :return:

getScatteringMatrixU(Ein, useTabulatedScatteringRadius=True, enableExtraCoulombPhase=False)[source]

Compute the scattering matrix

Note, unlike the getScatteringMatrixU() function in other resonance classes, the fact that different resonances are entirely different reactions means that the channelDicts have to have an additional layer of sorting that corresponds to the SLBW “level”.

setResonanceParametersByChannel(multipleSScheme='ENDF', useReichMooreApproximation=False, Ein=None, warnOnly=False)[source]

Reorganize member data into channels (relies heavily on groundwork in sortLandJ).

Note, unlike the getResonanceParametersByChannel() function in MLBW, RM or RML, the fact that different resonances are entirely different reactions means that the channelDicts have to have an additional layer of sorting that corresponds to the SLBW “level”.

class fudge.processing.resonances.reconstructResonances.URRPDFTable(*args, **kwds)[source]

Bases: collections.OrderedDict

Lightweight class for storage and IO of URR PDF’s

Test it as follows:

>>> import xs_pdf as xp
>>> import numpy as np
>>> u=xp.URRPDFTable()
>>> u.eBins=xp.equal_bins(2,0.,4.)
>>> u.xsBins=xp.equal_bins(3,1.,1.5)
>>> u['capture']=np.zeros(shape=(2,3))
>>> u['capture'][0,0]=1.
>>> u['capture'][1,2]=2.
>>> print(u)
    {"eBins": "[ 0.  2.  4.]", "urrPdf": {"capture": "[[ 1.  0.  0.]
[ 0. 0. 2.]]”}, “xsBins”: “[ 1. 1.16666667 1.33333333 1.5 ]”}
>>> u.save('junk.txt')
>>> v=xp.URRPDFTable()
>>> v.load('junk.txt')
>>> print(v)
    {"eBins": "[ 0.  2.  4.]", "urrPdf": {"capture": "[[ 1.  0.  0.]

[ 0. 0. 2.]]”}, “xsBins”: “[ 1. 1.16666667 1.33333333 1.5 ]”}

Emax
Emin
NE
NXS
XSmax
XSmin
eBinCenters
eBinWidths
load(fname)[source]
normalize()[source]

Normalize the binned PDFs using numpy trickery

pdf_at_energy(key, E)[source]
save(fname, verbose=False)[source]
xsBinCenters
xsBinWidths
class fudge.processing.resonances.reconstructResonances.URRcrossSection(reactionSuite, verbose=True, **kw)[source]

Bases: fudge.processing.resonances.reconstructResonances.resonanceReconstructionBaseClass

generateEnergyGrid(interpolateWidths=True)[source]

Get the energy grid for reconstructing unresolved resonances. Usually this is just the energy grid chosen by the evaluator for storing energy-dependent widths.

ENDF has been revised (with VIII) to interpolate on widths rather than reconstructed cross sections. For discussion see Red’s rant in D.E. Cullen “A Short History of ENDF/B Unresolved Resonance Parameters”, LLNL Report LLNL-TR-461199, ENDF Report ENDF-369 (2010)).

Parameters:interpolateWidths – if True, interpolate the average widths
Returns:a tuple containing the grid and the interpolate flag
getCrossSection(E, **kwargs)[source]
getFluctuationIntegrals(E, widths, DOF)[source]

From subroutine GNRL3 in RECENT. If possible, this will be replaced with more basic approach (rather than using lookup table)… not finding appropriate equations right now

Comments from GNRL3 sourcecode:

Calculate unresolved resonance fluctuation function
(original coding from AVERAGE4 by Mulki Bhat)
(new weighting scheme from MC^2-II)

This routine has been modified to calculate elastic, capture
and fission fluctuation functions all during one call (as
opposed to the original version which calculated each reaction
separately).

GNX, GGX, GFX and GXX are the widths for elastic, capture,
fission and competition. MUN, MUF and MUX are the number of
degrees of freedom for elastic, fission and competition (infinite
number of degrees assumed for capture). RN, RC and RF are the
calculated fluctuation integrals for elastic, capture and fission

The number of degrees of freedom for each distribution (elastic,
fission or competition) may be 1 to 4. If the number of degrees
of freedom for any distribution is less than 1 or more than 4
it will be treated as an infinite number of degrees of freedom
(which infers that the widths are not distributed, but are rather
all equal to the average value). This last case is simulated by
defining an additional 10 point quadrature in which the weight
for one point is 1.0 and the weight for all other points is zero.
for the one point of weight 1.0 the average width will be used.

:type E: Incident energy, may be scalar or numpy array
:type widths: dictionary of numpy arrays,
  containing average widths for each channel (for a specific L/J combination)
:type DOF: dictionary of degrees of freedom for specific L/J

:returns tuple(numpy.array) with fluctuation integrals for elastic, capture and fission
getLastResolvedResonanceEnergy(l=None, j=None)[source]

Get the last resonance energy from the resolved region, that will start all the ladders in the URR

Parameters:
  • l – orbital angular momentum of the resonance to find
  • j – total angular momentum of the resonance to find
Returns:

either the energy of the last resonance with requested (l,j) or None if it can’t be found

getLastResolvedResonanceRegion()[source]

Get the highest energy resonance ENDF region in the resolved resonance region (which may contain more than one ENDF region)

Returns:The ENDF resolved resonance region, or None if there is no ENDF region in the data
getTransmissionCoefficientsFromSumRule(skipFission=True)[source]
Use Moldauer’s sum rule to extract the transmission coefficients directly from the URR tables
P.A. Moldauer Phys. Rev. Lett. 19, 1047-1048 (1967)
Parameters:skipFission – flag to skip fission, what else?
Returns:a dictionary of results, sorted by channel
getURRPDF(resClass, nSamples=2000, style='goe', interpolateWidths=True, NE=100, ELowerBound=None, EUpperBound=None, NXS=200, XSLowerBound=None, XSUpperBound=None, verbose=True, restart=False, timing=False, memory=False)[source]
About restart files: they json dumps of the URR PDF, but before being renomalized. In other words, they are pure histograms
of number of points of the cross section curve in an energy/cross section voxel. This is a poor-mans path-length in the voxel.

At end of the run, the URR PDF (which, as I said, is still just a histogram) is normalized properly, turning it into a true PDF.

Parameters:
  • resClass – The class instance to use for reconstructing realizations of the cross section set.
  • nSamples
  • style – not the GNDS style, rather the style of resonance generation (see sampleRR)
  • interpolateWidths – interpolate the URR widths or not (just say “True”)
  • NE – number of energy bins in the PDF (# bin boundaries is NE+1)
  • ELowerBound
  • EUpperBound
  • NXS – number of cross section bins in the PDF (# bin boundaries is NXS+1)
  • XSLowerBound
  • XSUpperBound
  • verbose – enable verbose output
  • restart – enable restart capability
  • timing – enable printing out timings
  • memory – enable printing out memory usage of the resonance data structure
Returns:

getWidthsAndSpacings()[source]

URR tables give us the average resonance parameters, the average resonance spacing and the number of degrees of freedom, assuming that the widths are distributed by a chi^2 probability density function with the given number of degrees of freedom.

However, later we’ll need the level densities. Below we’ll construct them from the average resonance spacing for each L, J.

This function sets several member data as dicts (indexed by L & J) of data and interpolable functions:

- levelSpacings::

    ..math::
        D_{L,J}(E)

- levelDensities: derived from level spacings::

    ..math::

ho_{L,J}(E)=1/D_{L,J}(E)

  • averageWidths, assigned simple labels like ‘elastic’ and ‘capture’ for convenience:

    ..math::
        \Gamma_{L,J,c}
    
  • DOFs: degrees of freedom with same labels as averageWidths

  • reactionLabels: which reactions correspond to ‘elastic’, ‘capture’ and ‘fission’

return:None
rho(E)[source]

The dimensionless parameter rho:

..math::

ho = k(E) a

We always calculate channel radius for unresolved region according to the ENDF manual

param E:the incident neutron energy
return:rho
sampleRR(lastResonanceEnergies, lowerBound=None, upperBound=None, style='goe', verbose=True)[source]

Generate a sample of a resolved resonant set using the average URR parameters

Parameters:
  • lastResonanceEnergies – the energy of the last energy from the RRR
  • verbose – turn on verbose output
Returns:

the resonance set as a fudge.core.math.table

fudge.processing.resonances.reconstructResonances.blockwise(function)[source]
fudge.processing.resonances.reconstructResonances.getAllowedTotalSpins(L, S, useFactor2Trick=True)[source]

Returns a list of allowed J values from summing angular momenta L and S, where \vec{J}=\vec{L}+\vec{S} which implies |L-S| \leq J \leq L+S

The useFactor2Trick flag tells the routine whether we are summing real angular momenta or momenta * 2. If the useFactor2Trick flag is true, then momenta are really momenta*2, meaning they can be pure integers, even if the real momenta refer to 1/2-integer values (e.g. spin). The default is to useFactor2Trick because most C/C++/Fortran codes that compute angular momentum-stuff use the trick so they can use integer math. Also, it makes the use of the Python range() function possible.

fudge.processing.resonances.reconstructResonances.getR_S(E, Eres, captureWidth, widths, penetrabilities)[source]

Both versions of Reich_Moore formalisms (LRF=3 and 7) rely on building the complex matrix ‘R’. Here we break it up into the real component R and the imaginary component S, which represent symmetric and anti-symmetric scattering respectively

matrix elements R[i,j] and S[i,j] =

(summed over resonances) partialWidth[i]*partialWidth[j] * coefficient/(dE**2+captureWidth**2),

for the ith/jth open channel. For S, the coefficient is ‘dE’, for R ‘captureWidth’ and partialWidth[i] = widths[i] * penetrabilities[i]

( widths[i] is a row vector of resonance widths, and penetrabilities[i] is a column vector of the penetrability for each incident energy )

After computing R and S, invert using invertMatrices to find RI and SI such that (I+R+jS)*(I+RI+jSI) = I where I is the identity matrix, and j = sqrt(-1)

Incident energy dependence appears in both in E and the penetrabilities

fudge.processing.resonances.reconstructResonances.getResonanceReconstructionClass(formalism)[source]
fudge.processing.resonances.reconstructResonances.invertMatrices(R, S)[source]

find RI and SI such that (I+R+jS)*(I+RI+jSI) = I, where I is the identity matrix, and j = sqrt(-1)

for more info, see comments for subroutine FROBNS3 in recent.f

fudge.processing.resonances.reconstructResonances.reconstructAngularDistributions(reactionSuite, tolerance=None, verbose=False, energyUnit='eV')[source]

Reconstruct all pure two-body angular distributions from the resonance region. For MLBW and RM, that means ‘elastic’ channel only, but R-Matrix evaluations can support additional channels. SLBW cannot be used for angular distributions.

Returns a Python dict. They key is the reaction and the value is a reconstructed distributions.angular.XYs2d instance containing Legendre expansions for each incident energy.

fudge.processing.resonances.reconstructResonances.reconstructResonances(reactionSuite, tolerance=None, verbose=False, significantDigits=None, energyUnit='eV', disableUnresolvedWidthInterpolation=False)[source]

Reconstruct all resonance cross sections (resolved and unresolved) in reactionSuite, and add results together for full (resonance region) pointwise cross section.

Optional arguments: ‘tolerance’: fractional tolerance, used to refine interpolation grid.

e.g. if tolerance=0.001, points are added until lin-lin interpolation is good to 0.1% everywhere

‘verbose’: print status messages during reconstruction

‘significantDigits’: Controls how many digits can be used to represent the energy grid. For example,
if significantDigits=4 the resulting energy grid can contain 1.034, 1.035, 1.036 but not 1.0345. Using significantDigits=8 should allow data to be written back to ENDF-6 without loss of precision.

‘energyUnit’: unit to use for reconstruction. Resonance energies and widths will be converted to this unit.

‘disableUnresolvedWidthInterpolation’: set ‘True’ to reproduce old ENDF recommendation to interpolate
in cross sections rather than widths. That recommendation was reversed before ENDF-VIII release
fudge.processing.resonances.reconstructResonances.reconstructURRCrossSectionPDF(reactionSuite, tolerance=None, verbose=False, reconstructionScheme=None)[source]
class fudge.processing.resonances.reconstructResonances.resonanceReconstructionBaseClass(reactionSuite, energyUnit='eV', **kw)[source]

Bases: object

getCrossSection(E)[source]
k(energy)[source]
For an incident neutron, with energy in eV, the ENDF manual states :math:`k =

rac{sqrt{2m_n}}{hbar} rac{AWRI}{AWRI+1}sqrt{|E|}`

sqrt( 2 * neutronMass ) / hbar == 2.196807e-3 (eV*barn)**-1/2. Thus for energy in eV, k is in b**-1/2
penetrationFactor(L, rho)[source]
phi(L, rho)[source]
refineInterpolation(egrid, xsecs, tolerance=0.01, significantDigits=None)[source]

generateEnergyGrid may not give a fine enough grid to linearly interpolate to desired tolerance. My solution to that: for all consecutive points (x0,y0), (x1,y1) and (x2,y2) do a linear interpolation between (x0,y0) and (x2,y2). If the interpolation doesn’t agree with (x1,y1) within tolerance, subdivide up the region by adding two more calculated points. Iterate until interpolation agrees within tolerance.

This means that in the end we will have more points than required for given tolerance. The results can be thinned (thinning implemented in xData.XYs1d)

shiftFactor(L, rho)[source]

calculate shift factor used in SLBW and MLBW formalisms

fudge.processing.resonances.reconstructResonances.spins_equal(s1, s2)[source]

fudge.processing.resonances.setup module

Module contents