NormalizerGSF

class ompy.NormalizerGSF(*, normalizer_nld=None, nld=None, nld_model=None, alpha=None, gsf=None, path='saved_run/normalizers', regenerate=False, norm_pars=None)[source]

Bases: AbstractNormalizer

Normalize γSF to a given` <Γγ> (Gg)

Normalizes nld/gsf according to the transformation eq (3), Schiller2000:

gsf' = gsf * B * np.exp(alpha * Eg)
Variables:
  • gsf_in (Vector) – gsf to normalize

  • nld (Vector) – nld

  • normalizer_nld (NormalizerNLD) – NormalizerNLD to retriev parameters.

  • nld_model (Callable[..., Any]) – Model for nld above data of the from y = nld_model(E)

  • alpha (float) – tranformation parameter α

  • model_low (ExtrapolationModelLow) – Extrapolation model for the gsf at low energies (below data)

  • model_high (ExtrapolationModelHigh) – Extrapolation model for the gsf at high energies (above data)

  • norm_pars (NormalizationParameters) – Normalization parameters like experimental <Γγ>

  • method_Gg (str) – Method to use for the <Γγ> integral

  • path (Path) – The path save the results.

  • res (ResultsNormalized) – Results

Note

The prefered syntax is Normalizer(nld=…) If neither is given, the nld (and other parameters) can be explicity be set later by:

`normalizer.normalize(..., nld=...)`

or:

`normalizer.nld = ...`

In the later case you might have to send in a copy if it’s a mutable to ensure it is not changed.

Parameters:
  • normalizer_nld (Optional[NormalizerNLD], optional) – NormalizerNLD to retrieve parameters. If nld and/or nld_model are not set, they are taken from normalizer_nld.res in normalize.

  • nld (Optional[Vector], optional) – NLD. If not set it is taken from normalizer_nld.res in normalize.

  • nld_model (Optional[Callable[..., Any]], optional) – Model for nld above data of the from y = nld_model(E). If not set it is taken from normalizer_nld.res in normalize.

  • alpha (Optional[float], optional) – tranformation parameter α

  • gsf (Optional[Vector], optional) – gsf to normalize.

  • norm_pars (Optional[NormalizationParameters], optional) – Normalization parameters like experimental <Γγ>

Attributes Summary

LOG

Methods Summary

Gg_before_norm()

Compute <Γγ> before normalization

Gg_standard()

Compute normalization from <Γγ> (Gg) integral, the "standard" way

SpinSum(Ex, J[, Ltransfer])

Sum of spin distributions of the available states x 2

SpinSum_save_reload(SpinSum_args)

Reload SpinSum if computed with the same parameters before

extrapolate([gsf, E])

Extrapolate gsf using given models

fgsf(E, gsf, gsf_low, gsf_high)

Function composed of gsf and model, providing y = gsf(E)

fnld(E, nld, nld_model)

Function composed of nld and model, providing y = nld(E)

lnlike()

log likelihood

load([path])

Loads (pickeled) instance.

normalize(*[, gsf, normalizer_nld, alpha, ...])

Normalize gsf to a given <Γγ> (Gg).

plot([ax, add_label, add_figlegend, ...])

Plot the gsf and extrapolation normalization

plot_interactive()

Interactive plot to study the impact of different fit regions

save([path, overwrite])

Save (pickels) the instance

save_results_txt([path, nld, gsf, samples, ...])

Save results as txt

self_if_none(*args, **kwargs)

wrapper for lib.self_if_none

spin_dist(Ex, J)

Wrapper for ompy.SpinFunctions() curried with model and pars

Attributes Documentation

LOG = <Logger ompy.normalizer_gsf (WARNING)>

Methods Documentation

Gg_before_norm()[source]

Compute <Γγ> before normalization

Returns:

<Γγ> before normalization, in meV

Return type:

Gg (float)

Gg_standard()[source]

Compute normalization from <Γγ> (Gg) integral, the “standard” way

Equals “old” (normalization.f) version in the Spin sum get the normaliation, see eg. eq (21) and (26) in Larsen2011; but converted T to gsf.

Assumptions: s-wave (current implementation) and equal parity

# better format in shpinx
To derive the calculations, we start with
<Γγ(E,J,π)> = 1/(2π ρ(E,J,π)) ∑_XL ∑_Jf,πf ∫dEγ 2π Eγ³ gsf(Eγ)
                                                * ρ(E-Eγ, Jf, πf)

which leads to (eq 26)**

<Γγ> = 1 / (4π ρ(Sₙ, Jₜ± ½, πₜ)) ∫dEγ 2π Eγ³ gsf(Eγ) ρ(Sₙ-Eγ) spinsum
     = 1 / (2  ρ(Sₙ, Jₜ± ½, πₜ)) ∫dEγ    Eγ³ gsf(Eγ) ρ(Sₙ-Eγ) spinsum
    (= 1 / (ρ(Sₙ, Jₜ± ½, πₜ)) ∫dEγ Eγ³ gsf(Eγ) ρ(Sₙ-Eγ) spinsum(π) )
    (= D0 1 ∫dEγ Eγ³ gsf(Eγ) ρ(Sₙ-Eγ) spinsum(π) )

where the integral runs from 0 to Sₙ, and the spinsum selects the
available spins in dippole decays j=[-1,0,1]:
spinsum = ∑ⱼ g(Sₙ-Eγ, Jₜ± ½+j).

and <Γγ> is a shorthand for <Γγ(Sₙ, Jₜ± ½, πₜ)>.

When Jₜ≠0 where we see resonances from both states Jₜ± ½, but
often only <Γγ> is provided. We assume this is calculated as the
average over all width: (dropping Sₙ and πₜ in notation)
<Γγ> = ( N(Jₜ+½) <Γγ(Jₜ+ ½)> + N(Jₜ-½) <Γγ(Jₜ-½)> )
        / ( N(Jₜ+½) + N(Jₜ-½) ),
where N(J) is the number of levels with spin J. This is by ρ(J).
N(Jₜ+½) + N(Jₜ-½) is the same as ρ(Sₙ, Jₜ± ½, πₜ).

We can obtain ρ(Sₙ, Jₜ± ½, πₜ) from eq(19) via D0:
ρ(Sₙ, Jₜ± ½, πₜ) = ρ(Sₙ, Jₜ+ ½, πₜ) + ρ(Sₙ, Jₜ+ ½, πₜ) = 1/D₀
                 = ½ (ρ(Sₙ, Jₜ+ ½) + ρ(Sₙ, Jₜ+ ½))  [equi-parity]

Equi-parity means further that g(J,π) = 1/2 * g(J), for the calc.
above: spinsum(π) =  1/2 * spinsum.

For the Jₜ≠0 case, we can now rewrite:
<Γγ> = D0 * (ρ(Sₙ, Jₜ+ ½) * <Γγ(Jₜ+ ½)>
             + ρ(Sₙ, Jₜ- ½) * <Γγ(Jₜ- ½)>),
     = D0 * (  ∫dEγ Eγ³ gsf(Eγ) ρ(Sₙ-Eγ) spinsum(Jₜ- ½, π)
             + ∫dEγ Eγ³ gsf(Eγ) ρ(Sₙ-Eγ) spinsum(Jₜ+ ½, π) )
     = D0 ∫dEγ Eγ³ gsf(Eγ) ρ(Sₙ-Eγ) (  spinsum(Jₜ- ½, π)
                                     + spinsum(Jₜ+ ½, π))


** We set B to 1, so we get the <Γγ> before normalization. The
normalization B is then B = <Γγ>_exp/ <Γγ>_cal
Return type:

float

Returns:

Calculated Gg from gsf and nld

SpinSum(Ex, J, Ltransfer=0)[source]
Sum of spin distributions of the available states x 2

(assumes equiparity)

Note

∑_It ∑_Jr g(Jₜ- ½), where the first sum is over the available states in the residual ‘ nuclus assuming a angular momentum transfer l: It = Jₜ ± ½ ± l. The second sum assumes that we can reach states only b dipole radiaton, so the available final states are: Jr = It ± 1.

As g(J) is normalized as ∑_J g(J) = 1, and not ∑_J g(J, π) = ½ the sum calculated here is actually 2x the sum of avaiable states.

Parameters:
  • Ex (Union[float, np.ndarray]) – Excitation energy

  • J (float) – Target spin (in neutron capture reaction)

  • Ltransfer (float, optional) – Angular momentum transfer. Default is 0, corresponding to “swave” transfer.

Returns:

Sum of spin distributions. If Ex is

and array, this will be an array with the sum for each Ex.

Return type:

Union[float, np.ndarray]

SpinSum_save_reload(SpinSum_args)[source]

Reload SpinSum if computed with the same parameters before

Note

Note the most beautiful comparison, but we get a significant speedup in the multinest calculations

Parameters:

SpinSum_args – Arguments sent to SpinSum

Returns:

True, if args of SpinSum are the same as before,

and self.spincutModel and self.spincutPars are unchanged

Return type:

bool

extrapolate(gsf=None, E=[None, None])[source]

Extrapolate gsf using given models

Parameters:
  • gsf (Optional[Vector]) – If extrapolation is fit, it will be fit to this vector. Default is self._gsf.

  • E (optional) – extrapolation energies [Elow, Ehigh]

Return type:

Tuple[Vector, Vector]

Returns:

The extrapolated gSF-vectors on the form [low region, high region]

Raises:

ValueError if the models have any None variables.

static fgsf(E, gsf, gsf_low, gsf_high)[source]

Function composed of gsf and model, providing y = gsf(E)

Note

It will take the extrapolation where no exp. data is available.

Parameters:
  • E (ndarray) – Energies to evaluate

  • gsf (Vector) – GSF Vector for composition (mid part)

  • gsf_low (Vector) – GSF Vector for composition (lower part), usually an extrapolation

  • gsf_high (Vector) – GSF Vector for composition (higher part), usually an extrapolation

Returns:

Description

Return type:

ndarray

static fnld(E, nld, nld_model)[source]

Function composed of nld and model, providing y = nld(E)

Parameters:
  • E (ndarray) – Energies to evaluate

  • nld (Vector) – NLD Vector for composition

  • nld_model (Callable) – NLD model for composition

Returns:

Composite NLD

Return type:

ndarray

Raises:

ValueError – For low energies, the nld is not extrapolated. This may cause an error, in the calculation of Gg, when nld does not extend to 0. Please see https://github.com/oslocyclotronlab/ompy/issues/170 for more info.

lnlike()[source]

log likelihood

load(path=None)

Loads (pickeled) instance.

Such that it can be loaded if regenerate = False. Note that if any modifications of the __getstate__ method are present, these will effect what attributes are pickeled.

Parameters:

path (Union[str, Path, None]) – The path to the directoryto load file. If the value is None, ‘self.path’ will be used.

Raises:

FileNotFoundError – If file is not found

normalize(*, gsf=None, normalizer_nld=None, alpha=None, nld=None, nld_model=None, norm_pars=None, num=0)[source]

Normalize gsf to a given <Γγ> (Gg). Saves results to self.res.

Parameters:
  • normalizer_nld (Optional[NormalizerNLD], optional) – NormalizerNLD to retrieve parameters. If nld and/or nld_model are not set, they are taken from normalizer_nld.res in normalize.

  • nld (Optional[Vector], optional) – NLD. If not set it is taken from normalizer_nld.res in normalize.

  • nld_model (Optional[Callable[..., Any]], optional) – Model for nld above data of the from y = nld_model(E). If not set it is taken from normalizer_nld.res in normalize.

  • alpha (Optional[float], optional) – tranformation parameter α

  • gsf (Optional[Vector], optional) – gsf to normalize.

  • norm_pars (Optional[NormalizationParameters], optional) – Normalization parameters like experimental <Γγ>

  • num (Optional[int], optional) – Loop number, defaults to 0.

Return type:

None

plot(ax=None, *, add_label=True, add_figlegend=True, plot_fitregion=True, results=None, reset_color_cycle=True, **kwargs)[source]

Plot the gsf and extrapolation normalization

Parameters:
  • ax (optional) – The matplotlib axis to plot onto. Creates axis is not provided

  • add_label (bool, optional) – Defaults to True.

  • add_figlegend (bool, optional) – Defaults to True.

  • plot_fitregion (Optional[bool], optional) – Defaults to True.

  • results (ResultsNormalized, optional) – If provided, gsf and model are taken from here instead.

  • reset_color_cycle (Optional[bool], optional) – Defaults to True

  • **kwargs – Additional keyword arguments

Return type:

Tuple[Any, Any]

Returns:

fig, ax

plot_interactive()[source]

Interactive plot to study the impact of different fit regions

Note

  • This implementation is not the fastest, however helped to reduce

    the code quite a lot compared to slider_update

save(path=None, overwrite=True)

Save (pickels) the instance

Such that it can be loaded, and enabling the regenerate later.

Parameters:
  • path (Union[str, Path, None]) – The path to the save directory. If the value is None, ‘self.path’ will be used.

  • overwrite (bool) – Overwrite file if existent

save_results_txt(path=None, nld=None, gsf=None, samples=None, suffix=None)

Save results as txt

Uses a folder to save nld, gsf, and the samples (converted to an array)

Parameters:
  • path (Union[str, Path, None]) – The path to the save directory. If the value is None, ‘self.path’ will be used.

  • nld (Optional[Vector]) – (unnormalized) NLD to save

  • gsf (Optional[Vector]) – (unnormalized) gSF to save

  • samples (Optional[dict]) – samples to use for normalization

  • suffix (Optional[str]) – suffix to append to filename, eg. iteration number

self_if_none(*args, **kwargs)[source]

wrapper for lib.self_if_none

spin_dist(Ex, J)[source]

Wrapper for ompy.SpinFunctions() curried with model and pars