-
Notifications
You must be signed in to change notification settings - Fork 3
/
lnpi.py
436 lines (363 loc) · 12.5 KB
/
lnpi.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
r"""
Inverse temperature expansion of macrostate distribution (:mod:`~thermoextrap.lnpi`)
====================================================================================
This is used to extrapolate, in inverse temperature :math:`\beta = (k_{\rm B} T)^{-1}`, the macrostate distribution function :math:`\ln\Pi` from transition matrix Monte Carlo simulations.
See :ref:`examples/usage/basic/macrostate_dist_extrap:macrostate distribution extrapolation` for example usage.
"""
from __future__ import annotations
import warnings
from functools import lru_cache
from typing import Hashable, Sequence
import attrs
import cmomy
import numpy as np
import xarray as xr
# from attrs import converters as attc
from attrs import field
from attrs import validators as attv
from module_utilities import cached
from . import beta as beta_xpan
from .core._attrs_utils import cache_field, convert_dims_to_tuple
from .core.sputils import get_default_indexed, get_default_symbol
from .data import DataCallbackABC
from .docstrings import DOCFILLER_SHARED
from .models import Derivatives, ExtrapModel, SymFuncBase, SymSubs
docfiller_shared = DOCFILLER_SHARED.levels_to_top("cmomy", "xtrap", "beta").decorate
################################################################################
# lnPi correction stuff
################################################################################
class lnPi_func_central(SymFuncBase):
r"""
Special case of u_func_central.
For lnPi, have:
.. math::
\newcommand{\ave}[1]{\langle #1 \rangle}
(\ln \Pi)' = \frac{d \ln \Pi}{d \beta} = \mu N - \ave{u} + \ave{u - \mu N}_{\rm GC}
where :math:`\ave{}` and :math:`\ave{}_{\rm GC}` are the canonical and grand canonical (GC) ensemble averages.
We ignore the GC average term, as it does not depend on N. Note that this is not
necessarily the case for molecular systems.
So, the first derivative of this function is :func:`thermoextrap.beta.u_func_central`.
We consider only a correction of the form:
.. math::
(\ln\Pi)_{\text{energy}} = \ln\Pi - \beta \mu N = \ln Q - \ln \Xi
where :math:`Q\text{ and }\Xi` are the canonical and GC partition functions, respectively. thus,
.. math::
\begin{align*}
(\ln\Pi)_{\text{energy}}' &= - U \\
(\ln\Pi)_{\text{energy}}'' &= -U' \\
&\,\,\vdots
\end{align*}
"""
nargs = 1
u = get_default_symbol("u")
lnPi0 = get_default_symbol("lnPi0")
mudotN = get_default_symbol("mudotN")
@classmethod
def deriv_args(cls):
return [*beta_xpan.u_func_central.deriv_args(), cls.lnPi0, cls.mudotN]
def fdiff(self, argindex=1):
(beta,) = self.args
return self.mudotN - beta_xpan.u_func_central(beta)
@classmethod
def eval(cls, beta):
if beta is None:
return cls.lnPi0
return None
class lnPi_func_raw(SymFuncBase):
"""Raw moments version."""
nargs = 1
u = get_default_indexed("u")
lnPi0 = get_default_symbol("lnPi0")
mudotN = get_default_symbol("mudotN")
@classmethod
def deriv_args(cls):
return [*beta_xpan.u_func.deriv_args(), cls.lnPi0, cls.mudotN]
def fdiff(self, argindex=1):
(beta,) = self.args
return self.mudotN - beta_xpan.u_func(beta, 1)
@classmethod
def eval(cls, beta):
if beta is None:
return cls.lnPi0
return None
@lru_cache(5)
@docfiller_shared
def factory_derivatives(
name="lnPi",
n=None,
d=None,
xalpha=False,
central=False,
expand=True,
post_func=None,
):
"""
Expansion for ln(Pi/Pi_0) (ignore bad parts of stuff).
Parameters
----------
name : str, default='lnPi'
If name is `'lnPi'`, then get derivatives of lnPi.
Otherwise, get derivative object for general `X`.
{n_order}
{d_order}
{xalpha}
{central}
{expand}
{post_func}
Returns
-------
~thermoextrap.models.Derivatives
See Also
--------
thermoextrap.beta.factory_derivatives
"""
if name == "lnPi":
beta = get_default_symbol("beta")
func = lnPi_func_central(beta) if central else lnPi_func_raw(beta)
derivs = beta_xpan.SymDerivBeta(func=func, expand=expand, post_func=post_func)
exprs = SymSubs(
derivs, subs_all={derivs.beta: "None"}, expand=False, simplify=False
)
return Derivatives.from_sympy(exprs, args=derivs.args)
return beta_xpan.factory_derivatives(
name=name,
n=n,
d=d,
xalpha=xalpha,
central=central,
post_func=post_func,
expand=expand,
)
def _is_xr(name, x):
if not isinstance(x, xr.DataArray):
msg = f"{name} must be an xr.DataArray"
raise TypeError(msg)
return x
@attrs.define
class lnPiDataCallback(DataCallbackABC):
"""
Class to handle metadata callbacks for lnPi data.
Parameters
----------
lnPi0 : DataArray
Reference value of lnPi.
mu : DataArray
Value of chemical potential. Must have dimension ``dims_comp``.
dims_n : hashable or sequence of hashable
Dimension(s) for number of particle(s). That is, the dimensions of lnPi0 corresponding to particle number.
dims_comp : hashable
Dimension corresponding to components.
ncoords : DataArray, optional
Count of number of particles for given particle number (vector) and component.
Must have dimensions ``dims_comp`` and ``dims_n``.
allow_resample : bool, default=False
If True, allow simplified resampling of ``lnPi0`` data.
"""
# TODO(wpk): rename dims_comp to dim_comp.
#: lnPi data
lnPi0: xr.DataArray = field(validator=attv.instance_of(xr.DataArray))
#: Chemical potential
mu: xr.DataArray = field(validator=attv.instance_of(xr.DataArray))
#: Dimensions for particle number
dims_n: Hashable | Sequence[Hashable] = field(converter=convert_dims_to_tuple)
#: Dimensions for component
dims_comp: Hashable = field()
#: Particle number coordinates
ncoords: xr.DataArray = field(validator=attv.instance_of(xr.DataArray))
#: Flag to allow/disallow resampling of ``lnPi0``.
allow_resample: bool = field(default=False)
_cache: dict = cache_field()
# TODO(wpk): using dims_n, dims_comp naming because this is what is used in lnPi module
# def __init__(self, lnPi0, mu, dims_n, dims_comp, ncoords=None):
# if isinstance(dims_n, str):
# dims_n = [dims_n]
# self.dims_n = dims_n
# self.dims_comp = dims_comp
# self.lnPi0 = lnPi0
# self.mu = _is_xr("mu", mu)
# if ncoords is None:
# ncoords = self.get_ncoords()
# self.ncoords = _is_xr("ncoords", ncoords)
def check(self, data) -> None:
pass
@ncoords.default
def _set_default_ncoords(self):
# create grid
ncoords = np.meshgrid(
*tuple(self.lnPi0[x].values for x in self.dims_n), indexing="ij"
)
return xr.DataArray(np.array(ncoords), dims=(self.dims_comp, *self.dims_n))
@property
def lnPi0_ave(self):
if isinstance(self.lnPi0, xr.DataArray):
return self.lnPi0
# assume lnPi0 is an averaging object ala cmomy
return self.lnPi0.values
@cached.prop
def mudotN(self):
"""Dot product of `self.mu` and `self.ncoords`, reduces along `self.dims_comp`."""
return xr.dot(self.mu, self.ncoords, dims=self.dims_comp)
def resample(self, data, meta_kws=None, **kws):
"""Resample lnPi0 data."""
if not self.allow_resample:
msg = (
"Must set `self.allow_resample` to `True` to use resampling. "
"Resampling here is handled in an ad-hoc way, and should be "
"used with care."
)
raise ValueError(msg)
warnings.warn(
"'Correct' resampling of lnPi should be handled externally. "
"This resamples the average lnPi values. Instead, it is "
"recommended to resample based on collection matrices, and "
"construct lnPi values based on these.",
category=UserWarning,
stacklevel=2,
)
# wrap in xarray object:
dc = cmomy.xCentralMoments.from_vals(
self.lnPi0.expand_dims(dim="_new", axis=0),
axis="_new",
mom_dims="_mom",
mom=1,
)
# resample and reduce
dc, _ = dc.resample_and_reduce(**kws)
# return new object
return self.new_like(lnPi0=dc.values.sel(_mom=1))
def derivs_args(self, data, derivs_args):
return (*tuple(derivs_args), self.lnPi0_ave, self.mudotN)
# def factory_data_values(
# uv,
# order,
# lnPi,
# mu,
# dims_n,
# dims_comp,
# ncoords=None,
# central=False,
# skipna=False,
# rec_dim="rec",
# umom_dim="umom",
# xmom_dim="xmom",
# val_dims="val",
# rep_dim="rep",
# deriv_dim=None,
# chunk=None,
# compute=None,
# xv=None,
# x_is_u=True,
# **kws,
# ):
# """
# Factory function to produce a Data object
# Parameters
# ----------
# uv : array-like
# energy values. These are not averaged
# order : int
# highest umom_dim to calculate
# skipna : bool, default=False
# if True, skip `np.nan` values in creating averages.
# Can make some "big" calculations slow
# rec_dim, umom_dim, val_dim, rep_dim, deriv_dim : str
# names of record (i.e. time), umom_dim, value, replicate,
# and derivative (with respect to alpha)
# chunk : int or dict, optional
# If specified, perform chunking on resulting uv, xv arrays.
# If integer, chunk with {rec: chunk}
# otherwise, should be a mapping of form {dim_0: chunk_0, dim_1: chunk_1, ...}
# compute : bool, optional
# if compute is True, do compute averages greedily.
# if compute is False, and have done chunking, then defer calculation of averages (i.e., will be dask future objects).
# Default is to do greedy calculation
# constructor : 'val'
# kws : dict, optional
# extra arguments
# """
# if central:
# cls = DataValuesCentral
# else:
# cls = DataValues
# meta = lnPiDataCallback(
# lnPi0=lnPi, mu=mu, dims_n=dims_n, dims_comp=dims_comp, ncoords=ncoords
# )
# return cls.from_vals(
# uv=uv,
# xv=xv,
# order=order,
# skipna=skipna,
# rec_dim=rec_dim,
# umom_dim=umom_dim,
# val_dims=val_dims,
# rep_dim=rep_dim,
# deriv_dim=deriv_dim,
# chunk=chunk,
# compute=compute,
# x_is_u=x_is_u,
# meta=meta,
# **kws,
# )
# much more likely to have pre-aves here, but save that for the user
@docfiller_shared
def factory_extrapmodel_lnPi(
beta,
data,
*,
central=None,
order=None,
alpha_name="beta",
derivatives=None,
post_func=None,
derivatives_kws=None,
):
"""
Factory function to create Extrapolation model for beta expansion.
Parameters
----------
{beta}
data : object
Data object.
Should include :class:`lnPiDataCallback` object as well
order : int, optional
maximum order.
If not specified, default to `data.order + 1`
{central}
{post_func}
{alpha_name}
derivatives : :class:`thermoextrap.models.Derivatives`, optional
Derivatives object. If not passed, construct derivatives using :func:`thermoextrap.lnpi.factory_derivatives`.
derivates_kws : mapping, optional
Optional parameters to :func:`thermoextrap.lnpi.factory_derivatives`.
Returns
-------
extrapmodel : :class:`~thermoextrap.models.ExtrapModel`
See Also
--------
thermoextrap.lnpi.factory_derivatives
~thermoextrap.models.ExtrapModel
"""
if central is None:
central = data.central
if order is None:
order = data.order + 1
if central != data.central:
raise ValueError
if order > data.order + 1:
raise ValueError
if not data.x_is_u:
raise ValueError
if derivatives is None:
if derivatives_kws is None:
derivatives_kws = {}
derivatives = factory_derivatives(
name="lnPi", central=central, post_func=post_func, **derivatives_kws
)
return ExtrapModel(
alpha0=beta,
data=data,
derivatives=derivatives,
order=order,
# minus_log=mineus_log,
alpha_name=alpha_name,
)