Skip to content

Commit

Permalink
More new depend work
Browse files Browse the repository at this point in the history
  • Loading branch information
ceriottm committed Nov 19, 2023
1 parent 07eefc3 commit 7c052bb
Show file tree
Hide file tree
Showing 7 changed files with 76 additions and 66 deletions.
37 changes: 21 additions & 16 deletions ipi/engine/motion/constrained_dynamics.py
Original file line number Diff line number Diff line change
Expand Up @@ -79,17 +79,17 @@ def __init__(
"""

super(Dynamics, self).__init__(fixcom=fixcom, fixatoms=fixatoms)
dd(self).dt = depend_value(name="dt", value=timestep)
self._dt = depend_value(name="dt", value=timestep)

if thermostat is None:
self.thermostat = Thermostat()
else:
self.thermostat = thermostat

if nmts is None or len(nmts) == 0:
dd(self).nmts = depend_array(name="nmts", value=np.asarray([1], int))
self._nmts = depend_array(name="nmts", value=np.asarray([1], int))
else:
dd(self).nmts = depend_array(name="nmts", value=np.asarray(nmts, int))
self._nmts = depend_array(name="nmts", value=np.asarray(nmts, int))

if barostat is None:
self.barostat = Barostat()
Expand All @@ -109,7 +109,7 @@ def __init__(
)

# splitting mode for the integrators
dd(self).splitting = depend_value(name="splitting", value=splitting)
self._splitting = depend_value(name="splitting", value=splitting)

# The list of constraints coming from the input is an actual list of independent
# constraints, and should be treated as such
Expand Down Expand Up @@ -199,7 +199,7 @@ def __init__(self, constraint_list, dt=1.0):
self.constraint_list = constraint_list

# time step - will have to be linked to the dynamics time step
dd(self).dt = depend_value(name="dt", value=dt)
self._dt = depend_value(name="dt", value=dt)

def bind(self, beads):
if beads.nbeads > 1:
Expand Down Expand Up @@ -229,6 +229,9 @@ def proj_manifold(self, beads, stepsize=None, proj_p=True):
raise NotImplementedError()


inject_depend_properties(ConstrainedDynamics, ["dt", "nmts", "splitting"])


class ConstraintSolver(ConstraintSolverBase):
"""An implementation of a constraint solver that uses M-RATTLE to
impose constraints onto the momenta and a quasi-Newton method
Expand Down Expand Up @@ -362,21 +365,20 @@ def bind(self, motion):
super(ConstrainedIntegrator, self).bind(motion)

self.constraint_list = motion.constraint_list
dself = dd(self)
if motion.nsteps_geo is None:
dself.nsteps_geo = depend_value(name="nsteps_geo", value=1)
self._nsteps_geo = depend_value(name="nsteps_geo", value=1)
else:
dself.nsteps_geo = depend_value(name="nsteps_geo", value=motion.nsteps_geo)
print(self.nsteps_geo)
dself.qdt.add_dependency(dself.nsteps_geo)
self._nsteps_geo = depend_value(name="nsteps_geo", value=motion.nsteps_geo)

self._qdt.add_dependency(self._nsteps_geo)
if motion.nsteps_o is None:
dself.nsteps_o = depend_value(name="nsteps_o", value=1)
self._nsteps_o = depend_value(name="nsteps_o", value=1)
else:
dself.nsteps_o = depend_value(name="nsteps_o", value=motion.nsteps_o)
print(self.nsteps_o)
dself.tdt.add_dependency(dself.nsteps_o)
dd(self).csolver = motion.csolver
dpipe(dself.qdt, dd(self.csolver).dt)
self._nsteps_o = depend_value(name="nsteps_o", value=motion.nsteps_o)

self._tdt.add_dependency(self._nsteps_o)
self.csolver = motion.csolver
dpipe(self._qdt, self.csolver._dt)

def proj_cotangent(self):
self.csolver.proj_cotangent()
Expand All @@ -385,6 +387,9 @@ def proj_manifold(self):
self.csolver.proj_manifold()


inject_depend_properties(ConstrainedIntegrator, ["nsteps_geo", "nsteps_o"])


class NVEConstrainedIntegrator(ConstrainedIntegrator):
"""A propagator for constant-energy integration under a set of constraints.
Expand Down
68 changes: 36 additions & 32 deletions ipi/engine/motion/dynamics.py
Original file line number Diff line number Diff line change
Expand Up @@ -68,10 +68,9 @@ def __init__(
"""

super(Dynamics, self).__init__(fixcom=fixcom, fixatoms=fixatoms)
dself = dd(self) # noqa

# initialize time step. this is the master time step that covers a full time step
dd(self).dt = depend_value(name="dt", value=timestep)
self._dt = depend_value(name="dt", value=timestep)

if thermostat is None:
self.thermostat = Thermostat()
Expand All @@ -86,9 +85,9 @@ def __init__(
self.thermostat = thermostat

if nmts is None or len(nmts) == 0:
dd(self).nmts = depend_array(name="nmts", value=np.asarray([1], int))
self._nmts = depend_array(name="nmts", value=np.asarray([1], int))
else:
dd(self).nmts = depend_array(name="nmts", value=np.asarray(nmts, int))
self._nmts = depend_array(name="nmts", value=np.asarray(nmts, int))

if barostat is None:
self.barostat = Barostat()
Expand All @@ -113,7 +112,7 @@ def __init__(
self.integrator = DummyIntegrator()

# splitting mode for the integrators
dd(self).splitting = depend_value(name="splitting", value=splitting)
self._splitting = depend_value(name="splitting", value=splitting)

# constraints
self.fixcom = fixcom
Expand Down Expand Up @@ -161,27 +160,26 @@ def bind(self, ens, beads, nm, cell, bforce, prng, omaker):
)

# Strips off depend machinery for easier referencing.
dself = dd(self)
dthrm = dd(self.thermostat)
dbaro = dd(self.barostat)
dens = dd(self.ensemble)

# n times the temperature (for path integral partition function)
dself.ntemp = depend_value(
self._ntemp = depend_value(
name="ntemp", func=self.get_ntemp, dependencies=[dens.temp]
)

# fixed degrees of freedom count
fixdof = self.get_fixdof()

# first makes sure that the thermostat has the correct temperature and timestep, then proceeds with binding it.
dpipe(dself.ntemp, dthrm.temp)
dpipe(self._ntemp, dthrm.temp)

# depending on the kind, the thermostat might work in the normal mode or the bead representation.
self.thermostat.bind(beads=self.beads, nm=self.nm, prng=prng, fixdof=fixdof)

# first makes sure that the barostat has the correct stress and timestep, then proceeds with binding it.
dpipe(dself.ntemp, dbaro.temp)
dpipe(self._ntemp, dbaro.temp)
dpipe(dens.pext, dbaro.pext)
dpipe(dens.stressext, dbaro.stressext)
self.barostat.bind(
Expand Down Expand Up @@ -243,7 +241,10 @@ def step(self, step=None):
self.ensemble.time += self.dt # increments internal time


class DummyIntegrator(dobject):
inject_depend_properties(Dynamics, ["dt", "nmts", "splitting", "ntemp"])


class DummyIntegrator:
"""No-op integrator for (PI)MD"""

def __init__(self):
Expand Down Expand Up @@ -284,41 +285,38 @@ def bind(self, motion):
self.fixatoms = motion.fixatoms
self.enstype = motion.enstype

dself = dd(self)
dmotion = dd(motion)

# no need to dpipe these are really just references
dself.splitting = dmotion.splitting
dself.dt = dmotion.dt
dself.nmts = dmotion.nmts
self._splitting = motion._splitting
self._dt = motion._dt
self._nmts = motion._nmts

# total number of iteration in the inner-most MTS loop
dself.inmts = depend_value(name="inmts", func=lambda: np.prod(self.nmts))
dself.nmtslevels = depend_value(name="nmtslevels", func=lambda: len(self.nmts))
self._inmts = depend_value(name="inmts", func=lambda: np.prod(self.nmts))
self._nmtslevels = depend_value(name="nmtslevels", func=lambda: len(self.nmts))
# these are the time steps to be used for the different parts of the integrator
dself.qdt = depend_value(
self._qdt = depend_value(
name="qdt",
func=self.get_qdt,
dependencies=[dself.splitting, dself.dt, dself.inmts],
dependencies=[self._splitting, self._dt, self._inmts],
) # positions
dself.pdt = depend_array(
self._pdt = depend_array(
name="pdt",
func=self.get_pdt,
value=np.zeros(len(self.nmts)),
dependencies=[dself.splitting, dself.dt, dself.nmts],
dependencies=[self._splitting, self._dt, self._nmts],
) # momenta
dself.tdt = depend_value(
self._tdt = depend_value(
name="tdt",
func=self.get_tdt,
dependencies=[dself.splitting, dself.dt, dself.nmts],
dependencies=[self._splitting, self._dt, self._nmts],
) # thermostat

dpipe(dself.qdt, self.nm._dt)
dpipe(dself.dt, dd(self.barostat).dt)
dpipe(dself.qdt, dd(self.barostat).qdt)
dpipe(dself.pdt, dd(self.barostat).pdt)
dpipe(dself.tdt, dd(self.barostat).tdt)
dpipe(dself.tdt, dd(self.thermostat).dt)
dpipe(self._qdt, self.nm._dt)
dpipe(self._dt, dd(self.barostat).dt)
dpipe(self._qdt, dd(self.barostat).qdt)
dpipe(self._pdt, dd(self.barostat).pdt)
dpipe(self._tdt, dd(self.barostat).tdt)
dpipe(self._tdt, dd(self.thermostat).dt)

if motion.enstype == "sc" or motion.enstype == "scnpt":
# coefficients to get the (baseline) trotter to sc conversion
Expand Down Expand Up @@ -387,6 +385,12 @@ def pconstraints(self):
bp[self.fixatoms * 3 + 2] = 0.0


inject_depend_properties(
DummyIntegrator,
["splitting", "nmts", "dt", "inmts", "nmtslevels", "qdt", "pdt", "tdt"],
)


class NVEIntegrator(DummyIntegrator):

"""Integrator object for constant energy simulations.
Expand Down Expand Up @@ -680,8 +684,8 @@ def bind(self, mover):
"""

super(SCIntegrator, self).bind(mover)
self.ensemble.add_econs(dd(self.forces).potsc)
self.ensemble.add_xlpot(dd(self.forces).potsc)
self.ensemble.add_econs(self.forces._potsc)
self.ensemble.add_xlpot(self.forces._potsc)

def pstep(self, level=0):
"""Velocity Verlet monemtum propagator."""
Expand Down
12 changes: 6 additions & 6 deletions ipi/engine/motion/motion.py
Original file line number Diff line number Diff line change
Expand Up @@ -7,12 +7,10 @@

import numpy as np

from ipi.utils.depend import depend_value
from ipi.utils.depend import dd
from ipi.utils.depend import dobject
from ipi.utils.depend import depend_value, inject_depend_properties


class Motion(dobject):
class Motion:

"""Base motion calculation class.
Expand Down Expand Up @@ -43,8 +41,7 @@ def __init__(self, fixcom=False, fixatoms=None):
initial positions.
"""

dself = dd(self)
dself.dt = depend_value(name="dt", value=0.0)
self._dt = depend_value(name="dt", value=0.0)
self.fixcom = fixcom
if fixatoms is None:
self.fixatoms = np.zeros(0, int)
Expand Down Expand Up @@ -85,3 +82,6 @@ def step(self, step=None):
"""Dummy simulation time step which does nothing."""

pass


inject_depend_properties(Motion, "dt")
10 changes: 6 additions & 4 deletions ipi/engine/motion/multi.py
Original file line number Diff line number Diff line change
Expand Up @@ -5,7 +5,7 @@
# i-PI Copyright (C) 2014-2015 i-PI developers
# See the "licenses" directory for full license information.

from ipi.utils.depend import depend_value, dd
from ipi.utils.depend import depend_value, inject_depend_properties
from ipi.engine.motion import Motion


Expand All @@ -21,11 +21,10 @@ def __init__(self, motionlist=None):
motion will be constrained or not. Defaults to False.
"""

dself = dd(self)
dself.dt = depend_value(name="dt", func=self.get_totdt)
dself._dt = depend_value(name="dt", func=self.get_totdt)

Check warning on line 24 in ipi/engine/motion/multi.py

View workflow job for this annotation

GitHub Actions / lint

F821 undefined name 'dself'
self.mlist = motionlist
for m in self.mlist:
dd(m).dt.add_dependant(dself.dt)
m._dt.add_dependant(self._dt)
a = ( # noqa
self.dt
) # DON'T ASK WHY BUT IF YOU DON'T DO THAT WEAKREFS TO SELF.DT WILL BE INVALIDATED
Expand Down Expand Up @@ -70,3 +69,6 @@ def bind(self, ens, beads, nm, cell, bforce, prng, omaker):

for m in self.mlist:
m.bind(ens, beads, nm, cell, bforce, prng, omaker)


inject_depend_properties(MultiMotion, "dt")
8 changes: 4 additions & 4 deletions ipi/engine/motion/planetary.py
Original file line number Diff line number Diff line change
Expand Up @@ -70,12 +70,9 @@ def __init__(
self.nbeads = nbeads
self.screen = screen

dself = dd(self)

# the planetary step just computes constrained-centroid properties so it
# should not advance the timer
dself.dt = depend_value(name="dt", value=0.0)
# dset(self, "dt", depend_value(name="dt", value = 0.0) )
self._dt = depend_value(name="dt", value=0.0)
self.fixatoms = np.asarray([])
self.fixcom = True

Expand Down Expand Up @@ -279,3 +276,6 @@ def step(self, step=None):
% (self.tmc / self.neval, self.tmtx / self.neval, self.tsave / self.neval),
verbosity.high,
)


inject_depend_properties(Planetary, "dt")
4 changes: 2 additions & 2 deletions ipi/engine/normalmodes.py
Original file line number Diff line number Diff line change
Expand Up @@ -167,7 +167,7 @@ def bind(self, ensemble, motion, beads=None, forces=None):

# stores a reference to the bound beads and ensemble objects
self.ensemble = ensemble
dpipe(dd(motion).dt, self._dt)
dpipe(motion._dt, self._dt)

# sets up what's necessary to perform nm transformation.
if self.nbeads == 1: # classical trajectory! don't waste time doing anything!
Expand Down Expand Up @@ -302,7 +302,7 @@ def bind(self, ensemble, motion, beads=None, forces=None):
)

self._dt = depend_value(name="dt", value=1.0)
dpipe(dd(self.motion).dt, self._dt)
dpipe(self.motion._dt, self._dt)
self._prop_pq = depend_array(
name="prop_pq",
value=np.zeros((self.beads.nbeads, 2, 2)),
Expand Down
3 changes: 1 addition & 2 deletions ipi/engine/smotion/metad.py
Original file line number Diff line number Diff line change
Expand Up @@ -8,7 +8,6 @@
# See the "licenses" directory for full license information.

from ipi.engine.smotion import Smotion
from ipi.utils.depend import *


__all__ = ["MetaDyn"]
Expand Down Expand Up @@ -90,7 +89,7 @@ def step(self, step=None):
for fc in s.ensemble.bias.mforces:
if fc.ffield == k:
for fb in fc._forces:
dd(fb).ufvx.taint()
fb._ufvx.taint()
meta_pot_after = s.ensemble.bias.pot
# updates the conserved quantity with the change in bias so that
# we remove the shift due to added hills
Expand Down

0 comments on commit 7c052bb

Please sign in to comment.