forcefields

Details

All pre-defined forcefield classes for use in flowerMD.

class flowermd.library.forcefields.BeadSpring(r_cut, beads, bonds=None, angles=None, dihedrals=None, exclusions=['bond', '1-3'])

Bead-spring forcefield class.

Given a dictionary of bead types, this class creates a list hoomd.md.force.Force objects to capture bonded and non-bonded interactions between the beads. For non-bonded interactions, a Lennard-Jones potential is used. For bonds and angles, a harmonic potential is used. For dihedrals, a periodic potential is used.

Parameters:
  • r_cut (float, required) – The cutoff radius for the LJ potential.

  • beads (dict, required) – A dictionary of bead types. Each bead type should be a dictionary with the keys “epsilon” and “sigma” that correspond to the LJ parameters.

  • bonds (dict, default None) – A dictionary of bond types separated by a dash. Each bond type should be a dictionary with the keys “r0” and “k” that correspond to the harmonic bond parameters.

  • angles (dict, default None) – A dictionary of angle types separated by a dash. Each angle type should be a dictionary with the keys “t0” and “k” that correspond to the harmonic angle parameters.

  • dihedrals (dict, default None) – A dictionary of dihedral types separated by a dash. Each dihedral type should be a dictionary with the keys “phi0”, “k”, “d”, and “n” that correspond to the periodic dihedral parameters.

  • exclusions (list, default ["bond", "1-3"]) – A list of exclusions to use in the neighbor list. The default is to exclude bonded and 1-3 interactions.

Examples

For a simple bead-spring model with two bead types A and B, the following code can be used:

ff = BeadSpring(r_cut=2.5,
        beads={"A": dict(epsilon=1.0, sigma=1.0),
               "B": dict(epsilon=2.0, sigma=2.0)},
        bonds={"A-A": dict(r0=1.1, k=300), "A-B": dict(r0=1.1, k=300)},
        angles={"A-A-A": dict(t0=2.0, k=200),
                "A-B-A": dict(t0=2.0, k=200)},
        dihedrals={"A-A-A-A": dict(phi0=0.0, k=100, d=-1, n=1)})
class flowermd.library.forcefields.EllipsoidForcefield(epsilon, lpar, lperp, r_cut, bond_k, bond_r0, angle_k=None, angle_theta0=None)

A forcefield for modeling anisotropic bead polymers.

Notes

This is designed to be used with flowermd.library.polymers.EllipsoidChain and uses ghost particles of type “A” and “B” for intra-molecular interactions of bonds and two-body angles. Ellipsoid centers (type “R”) are used in inter-molecular pair interations.

The set of interactions are: 1. hoomd.md.bond.Harmonic: Models ellipsoid bonds as tip-to-tip bonds 2. hoomd.md.angle.Harmonic: Models angles of two neighboring ellipsoids. 3. hoomd.md.pair.aniso.GayBerne” Model pair interactions between beads.

Parameters:
  • epsilon (float, required) – energy

  • lpar (float, required) – Semi-axis length of the ellipsoid along the major axis.

  • lperp (float, required) – Semi-axis length of the ellipsoid along the minor axis.

  • r_cut (float, required) – Cut off radius for pair interactions

  • bond_k (float, required) – Spring constant in harmonic bond

  • bond_r0 (float, required) – Equilibrium tip-to-tip bond length.

  • angle_k (float, required) – Spring constant in harmonic angle.

  • angle_theta0 (float, required) – Equilibrium angle between 2 consecutive beads.

class flowermd.library.forcefields.FF_from_file(*args: Any, **kwargs: Any)

Forcefield class for loading a forcefield from an XML file.

class flowermd.library.forcefields.GAFF(*args: Any, **kwargs: Any)

General Amber forcefield class.

class flowermd.library.forcefields.KremerGrestBeadSpring(bond_k, bond_max, radial_shift=0, sigma=1.0, epsilon=1.0, bead_name='A')

Kremer-Grest Bead-Spring polymer coarse-grain model.

Parameters:
  • bond_k (float, required) – Spring constant in the FENE-WCA bond potential.

  • bond_max (float, required) – Maximum bond length in the FENE-WCA bond potential.

  • delta (float, optional, default 0.0) – The radial shift used in the FENE-WCA bond potential.

  • sigma (float, optional, default 1.0) – Length scale in the 12-6 Lennard-Jones pair force.

  • epsilon (float, optional, default 1.0) – Energy scale in the 12-6 Lennard-Jones pair force.

  • bead_name (str, optional, default "A") – Particle names in the bead-spring system.

Notes

Use this forcefield class with flowermd.library.polymers.BeadSpring.

This forcefield class returns two types of interactions:

  1. 12-6 LJ pair potential with a cutoff of \(2^{(1/6)}\sigma\).

  2. Bond potential that includes a FENE spring and a WCA repulsive term.

The sigma and epsilon parameters are used both for the repulsive LJ potential and the WCA part of the bond potential.

class flowermd.library.forcefields.OPLS_AA(*args: Any, **kwargs: Any)

OPLS All Atom forcefield class.

class flowermd.library.forcefields.OPLS_AA_BENZENE(*args: Any, **kwargs: Any)

OPLS All Atom for benzene molecule forcefield class.

class flowermd.library.forcefields.OPLS_AA_DIMETHYLETHER(*args: Any, **kwargs: Any)

OPLS All Atom for dimethyl ether molecule forcefield class.

class flowermd.library.forcefields.OPLS_AA_PPS(*args: Any, **kwargs: Any)

OPLS All Atom for PPS molecule forcefield class.

class flowermd.library.forcefields.TableForcefield(pairs=None, bonds=None, angles=None, dihedrals=None, r_min=None, r_cut=None, exclusions=['bond', '1-3'], nlist_buffer=0.4)

Create a set of hoomd table potentials.

This class provides an interface for creating hoomd table potentials either from arrays of energy and forces, or from files storing the tabulated energy and forces.

In HOOMD-Blue, table potentials are available for:

  • Pairs: hoomd.md.pair.Table

  • Bonds: hoomd.md.bond.Table

  • Angles: hoomd.md.angle.Table

  • Dihedrals: hoomd.md.dihedral.Table

Notes

HOOMD table potentials are initialized using arrays of energy and forces. It may be most convenient to store tabulated data in files, in that case use the from_files method.

Parameters:
  • pairs (dict, optional, default None) –

  • bonds (dict, optional, default None) –

  • angles (dict, optional, default None) –

  • dihedrals (dict, optional, default None) –

  • r_min (float, optional, default None) – Sets the r_min value for hoomd.md.pair.Table parameters.

  • r_max (float, optional, default None) – Sets the r cutoff value for hoomd.md.pair.Table parameters.

  • exclusions (list of str, optional, default ["bond", "1-3"]) –

    Sets exclusions for hoomd.md.pair.Table neighbor list.

    See documentation for hoomd.md.nlist # noqa: E501

classmethod from_files(pairs=None, bonds=None, angles=None, dihedrals=None, exclusions=['bond', '1-3'], nlist_buffer=0.4, **kwargs)

Create a table forefield from files.

Parameters:
  • pairs (dict, optional, default None) – Dictionary with keys of pair type and keys of file path

  • bonds (dict, optional, default None) – Dictionary with keys of bond type and keys of file path

  • angles (dict, optional, default None) – Dictionary with keys of angle type and keys of file path

  • dihedrals (dict, optional, default None) – Dictionary with keys of dihedral type and keys of file path

  • **kwargs (keyword arguments) – Key word arguments passed to numpy.genfromtxt or numpy.load

Notes

The parameters must use a {“type”: “file_path”} mapping.

Following HOOMD conventions, pair types must be given as a tuple of particles types while bonds, angles and dihedrals are given as a str of particle types separated by dashes.

Example

table_forcefield = TableForcefield.from_files(
    pairs = {
        ("A", "A"): "A_pairs.txt"
        ("B", "B"): "B_pairs.txt"
        ("A", "B"): "AB_pairs.txt"
    },
    bonds = {"A-A": "A_bonds.txt", "B-B": "B_bonds.txt"},
    angles = {"A-A-A": "A_angles.txt", "B-B-B": "B_angles.txt"},
)

Warning

It is assumed that the structure of the files are:
  • Column 1: Independent variable (e.g. distance, length, angle)

  • Column 2: Energy

  • Column 3: Force