Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

get rid of npdim #457

Merged
merged 2 commits into from
Nov 6, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
2 changes: 2 additions & 0 deletions CHANGELOG.md
Original file line number Diff line number Diff line change
Expand Up @@ -7,7 +7,9 @@ The format is based on [Keep a Changelog](https://keepachangelog.com/en/1.0.0/),

[Unreleased]
### Added

### Changed
- Get rid of npdim option that at some point may have allowed the prior transformation to return higher dimensional vector than the inputs. Note that due to this change, restoring the checkpoint from previous version of the dynesty won't be possible) (issues #456, #457) (original issue reported by @MichaelDAlbrow, fixed by @segasai )
### Fixed


Expand Down
44 changes: 22 additions & 22 deletions py/dynesty/dynamicsampler.py
Original file line number Diff line number Diff line change
Expand Up @@ -345,7 +345,7 @@ def _initialize_live_points(live_points,
loglikelihood,
M,
nlive=None,
npdim=None,
ndim=None,
rstate=None,
blob=False,
use_pool_ptform=None):
Expand All @@ -372,7 +372,7 @@ def _initialize_live_points(live_points,
nlive: int
Number of live-points

npdim: int
ndim: int
Number of dimensions

rstate: :class: numpy.random.RandomGenerator
Expand Down Expand Up @@ -403,13 +403,13 @@ def _initialize_live_points(live_points,
# sampling from the unit cube.
n_attempts = 1000

min_npoints = min(nlive, max(npdim + 1, min(nlive - 20, 100)))
min_npoints = min(nlive, max(ndim + 1, min(nlive - 20, 100)))
# the minimum number points we want with finite logl
# we want want at least npdim+1, because we want
# we want want at least ndim+1, because we want
# to be able to constraint the ellipsoid
# Note that if nlive <npdim+ 1 this doesn't really make sense
# Note that if nlive <ndim+ 1 this doesn't really make sense
# but we should have warned the user earlier, so they are on their own
# And the reason we have max(npdim+1, X ) is that we'd like to get at
# And the reason we have max(ndim+1, X ) is that we'd like to get at
# least X points as otherwise the poisson estimate of the volume will
# be too large.
# The reason why X is min(nlive-20, 100) is that we want at least 100
Expand All @@ -419,8 +419,8 @@ def _initialize_live_points(live_points,
# integrals and batch sampling in plateau edge tests
# The formula probably should be simplified

live_u = np.zeros((nlive, npdim))
live_v = np.zeros((nlive, npdim))
live_u = np.zeros((nlive, ndim))
live_v = np.zeros((nlive, ndim))
live_logl = np.zeros(nlive)
ngoods = 0 # counter for how many finite logl we have found
live_blobs = []
Expand All @@ -429,7 +429,7 @@ def _initialize_live_points(live_points,
iattempt += 1

# simulate nlive points by uniform sampling
cur_live_u = rstate.random(size=(nlive, npdim))
cur_live_u = rstate.random(size=(nlive, ndim))
if use_pool_ptform:
cur_live_v = M(prior_transform, np.asarray(cur_live_u))
else:
Expand Down Expand Up @@ -611,7 +611,7 @@ def _configure_batch_sampler(main_sampler,
batch_sampler = _SAMPLERS[main_sampler.bounding](
main_sampler.loglikelihood,
main_sampler.prior_transform,
main_sampler.npdim,
main_sampler.ndim,
main_sampler.live_init, # this is not used at all
# as we replace the starting points
main_sampler.method,
Expand Down Expand Up @@ -658,7 +658,7 @@ def _configure_batch_sampler(main_sampler,
main_sampler.loglikelihood,
main_sampler.M,
nlive=nlive_new,
npdim=main_sampler.npdim,
ndim=main_sampler.ndim,
rstate=main_sampler.rstate,
blob=main_sampler.blob,
use_pool_ptform=main_sampler.use_pool_ptform)
Expand Down Expand Up @@ -776,7 +776,7 @@ def _configure_batch_sampler(main_sampler,
# Trigger an update of the internal bounding distribution based
# on the "new" set of live points.

live_u = np.empty((nlive_new, main_sampler.npdim))
live_u = np.empty((nlive_new, main_sampler.ndim))
live_v = np.empty((nlive_new, saved_v.shape[1]))
live_logl = np.empty(nlive_new)
live_bound = np.zeros(nlive_new, dtype=int)
Expand Down Expand Up @@ -867,7 +867,7 @@ class DynamicSampler:
Function transforming a sample from the a unit cube to the parameter
space of interest according to the prior.

npdim : int, optional
ndim : int, optional
Number of parameters accepted by `prior_transform`.

bound : {`'none'`, `'single'`, `'multi'`, `'balls'`, `'cubes'`}, optional
Expand Down Expand Up @@ -913,14 +913,14 @@ class DynamicSampler:
A dictionary of additional parameters (described below).
"""

def __init__(self, loglikelihood, prior_transform, npdim, bound, method,
def __init__(self, loglikelihood, prior_transform, ndim, bound, method,
update_interval_ratio, first_update, rstate, queue_size, pool,
use_pool, ncdim, nlive0, kwargs):

# distributions
self.loglikelihood = loglikelihood
self.prior_transform = prior_transform
self.npdim = npdim
self.ndim = ndim
self.ncdim = ncdim
self.blob = kwargs.get('blob') or False
# bounding/sampling
Expand Down Expand Up @@ -1205,7 +1205,7 @@ def sample_initial(self,
Contains `live_u`, the coordinates on the unit cube, `live_v`, the
transformed variables, and `live_logl`, the associated
loglikelihoods. By default, if these are not provided the initial
set of live points will be drawn from the unit `npdim`-cube.
set of live points will be drawn from the unit `ndim`-cube.
**WARNING: It is crucial that the initial set of live points have
been sampled from the prior. Failure to provide a set of valid
live points will lead to incorrect results.**
Expand All @@ -1216,7 +1216,7 @@ def sample_initial(self,
Index of the live point with the worst likelihood. This is our
new dead point sample.

ustar : `~numpy.ndarray` with shape (npdim,)
ustar : `~numpy.ndarray` with shape (ndim,)
Position of the sample.

vstar : `~numpy.ndarray` with shape (ndim,)
Expand Down Expand Up @@ -1292,7 +1292,7 @@ def sample_initial(self,
self.loglikelihood,
self.M,
nlive=nlive,
npdim=self.npdim,
ndim=self.ndim,
rstate=self.rstate,
blob=self.blob,
use_pool_ptform=self.use_pool_ptform)
Expand All @@ -1317,7 +1317,7 @@ def sample_initial(self,
first_update = self.first_update
self.sampler = _SAMPLERS[bounding](self.loglikelihood,
self.prior_transform,
self.npdim,
self.ndim,
self.live_init,
self.method,
update_interval,
Expand Down Expand Up @@ -1502,7 +1502,7 @@ def sample_batch(self,
new dead point sample. **Negative values indicate the index
of a new live point generated when initializing a new batch.**

ustar : `~numpy.ndarray` with shape (npdim,)
ustar : `~numpy.ndarray` with shape (ndim,)
Position of the sample.

vstar : `~numpy.ndarray` with shape (ndim,)
Expand Down Expand Up @@ -1972,7 +1972,7 @@ def run_nested(self,
Contains `live_u`, the coordinates on the unit cube, `live_v`, the
transformed variables, and `live_logl`, the associated
loglikelihoods. By default, if these are not provided the initial
set of live points will be drawn from the unit `npdim`-cube.
set of live points will be drawn from the unit `ndim`-cube.
**WARNING: It is crucial that the initial set of live points have
been sampled from the prior. Failure to provide a set of valid
live points will lead to incorrect results.**
Expand Down Expand Up @@ -2028,7 +2028,7 @@ def run_nested(self,
# The reason to scale with square of number of
# dimensions is because the number coefficients
# defining covariance is roughly 0.5 * N^2
n_effective = max(self.npdim * self.npdim, 10000)
n_effective = max(self.ndim * self.ndim, 10000)

stop_kwargs['target_n_effective'] = n_effective
nlive_init = nlive_init or self.nlive0
Expand Down
55 changes: 30 additions & 25 deletions py/dynesty/dynesty.py
Original file line number Diff line number Diff line change
Expand Up @@ -329,13 +329,6 @@ def prior_transform(u):
efficiency in percent (`'min_eff'`). Defaults are `2 * nlive` and
`10.`, respectively.

npdim : int, optional
Number of parameters accepted by `prior_transform`. This might
differ from `ndim` in the case where a parameter of loglikelihood
is dependent upon multiple independently distributed parameters,
some of which may
be nuisance parameters.

rstate : `~numpy.random.Generator`, optional
`~numpy.random.Generator` instance. If not given, the
global random state of the `~numpy.random` module will be used.
Expand Down Expand Up @@ -370,7 +363,7 @@ def prior_transform(u):
transformed variables, and `live_logl`, the associated
loglikelihoods.
By default, if these are not provided the initial set of live
points will be drawn uniformly from the unit `npdim`-cube.
points will be drawn uniformly from the unit `ndim`-cube.
**WARNING: It is crucial that the initial set of live points have
been sampled from the prior. Failure to provide a set of valid
live points
Expand Down Expand Up @@ -468,13 +461,16 @@ def prior_transform(u):
will be sampled using the sampling method, the remaining
dimensions will
just sample uniformly from the prior distribution.
If this is `None` (default), this will default to npdim.
If this is `None` (default), this will default to ndim.

blob: bool, optional
The default value is False. If it is true, then the log-likelihood
should return the tuple of logl and a numpy-array "blob" that will
stored as part of the chain. That blob can contain auxiliary
information computed inside the likelihood function.

npdim : int
This option is deprecated and should not be used
"""

static_docstring = f"""
Expand Down Expand Up @@ -548,10 +544,15 @@ def __new__(cls,
history_filename=None):

# Prior dimensions.
if npdim is None:
npdim = ndim
if ncdim is None:
ncdim = npdim
if npdim is not None:
if npdim != ndim:
raise ValueError('''npdim functionality is not functioning
and is deprecated ''')
else:
warnings.warn(
"""the npdim keyword/functionality is deprecated as not
functioning and will be removed in further releases""", DeprecationWarning)
ncdim = ncdim or ndim

# Bounding method.
if bound not in _SAMPLERS:
Expand All @@ -563,7 +564,7 @@ def __new__(cls,

walks, slices = _get_walks_slices(walks, slices, sample, ndim)

if ncdim != npdim and sample in ['slice', 'hslice', 'rslice']:
if ncdim != ndim and sample in ['slice', 'hslice', 'rslice']:
raise ValueError('ncdim unsupported for slice sampling')

# Custom sampling function.
Expand All @@ -585,7 +586,7 @@ def __new__(cls,
warnings.warn(
"Beware! Having `nlive <= 2 * ndim` is extremely risky!")

nonbounded = get_nonbounded(npdim, periodic, reflective)
nonbounded = get_nonbounded(ndim, periodic, reflective)
kwargs['nonbounded'] = nonbounded
kwargs['periodic'] = periodic
kwargs['reflective'] = reflective
Expand Down Expand Up @@ -680,7 +681,7 @@ def __new__(cls,
loglike,
M,
nlive=nlive,
npdim=npdim,
ndim=ndim,
rstate=rstate,
blob=blob,
use_pool_ptform=use_pool.get('prior_transform', True))
Expand All @@ -689,7 +690,7 @@ def __new__(cls,
sampler = super().__new__(_SAMPLERS[bound])
sampler.__init__(loglike,
ptform,
npdim,
ndim,
live_points,
sample,
update_interval,
Expand Down Expand Up @@ -754,11 +755,15 @@ def __init__(self,
history_filename=None):

# Prior dimensions.
if npdim is None:
npdim = ndim
if ncdim is None:
ncdim = npdim

if npdim is not None:
if npdim != ndim:
raise ValueError('''npdim functionality is not functioning
and is deprecated ''')
else:
warnings.warn(
"""the npdim keyword/functionality is deprecated as not
functioning and will be removed in further releases""", DeprecationWarning)
ncdim = ncdim or ndim
nlive = nlive or 500

# Bounding method.
Expand All @@ -771,7 +776,7 @@ def __init__(self,

walks, slices = _get_walks_slices(walks, slices, sample, ndim)

if ncdim != npdim and sample in ['slice', 'hslice', 'rslice']:
if ncdim != ndim and sample in ['slice', 'hslice', 'rslice']:
raise ValueError('ncdim unsupported for slice sampling')

update_interval_ratio = _get_update_interval_ratio(
Expand All @@ -791,7 +796,7 @@ def __init__(self,
# Citation generator.
kwargs['cite'] = _get_citations('dynamic', bound, sample)

nonbounded = get_nonbounded(npdim, periodic, reflective)
nonbounded = get_nonbounded(ndim, periodic, reflective)
kwargs['nonbounded'] = nonbounded
kwargs['periodic'] = periodic
kwargs['reflective'] = reflective
Expand Down Expand Up @@ -877,7 +882,7 @@ def __init__(self,
kwargs['compute_jac'] = compute_jac

# Initialize our nested sampler.
super().__init__(loglike, ptform, npdim, bound, sample,
super().__init__(loglike, ptform, ndim, bound, sample,
update_interval_ratio, first_update, rstate,
queue_size, pool, use_pool, ncdim, nlive, kwargs)

Expand Down