smc_lammps.generate package

Subpackages

Submodules

smc_lammps.generate.default_parameters module

class smc_lammps.generate.default_parameters.Parameters(N: int = 501, cycles: int | None = 2, max_steps: int | None = None, RNA_polymerase_size: float = 5.0, spaced_beads_interval: int | None = None, spaced_beads_size: float = 5.0, spaced_beads_full_dna: bool = False, spaced_beads_smc_clearance: float = 5.0, spaced_beads_custom_stiffness: float = 1.0, add_side_site: bool = False, site_cycle_period: int = 0, site_toggle_delay: int = 0, site_cycle_when: str = 'apo')

Bases: object

Class that stores all simulation parameters defined by the user.

N: int = 501

Number of DNA beads

RNA_polymerase_size: float = 5.0

radius of RNA polymerase (nm)

RNA_polymerase_type = 1
T = 300.0

Simulation temperature (K)

__init__(N: int = 501, cycles: int | None = 2, max_steps: int | None = None, RNA_polymerase_size: float = 5.0, spaced_beads_interval: int | None = None, spaced_beads_size: float = 5.0, spaced_beads_full_dna: bool = False, spaced_beads_smc_clearance: float = 5.0, spaced_beads_custom_stiffness: float = 1.0, add_side_site: bool = False, site_cycle_period: int = 0, site_toggle_delay: int = 0, site_cycle_when: str = 'apo') None
add_RNA_polymerase = True

adds 10 nm bead at DNA-tether site

add_side_site: bool = False

Add a binding site on the lower SMC arm to act as the cycling site. If enabled, the lower site operates normally.

add_stopper_bead = False

add a bead that prevents the SMC from slipping off of the wrong end of the DNA

arm_length = 50.0

Length of each coiled-coil arm (nm)

arms_angle_ATP = 130.0

Opening angle of arms in ATP-bound state (degrees)

arms_stiffness = 100.0

Bending stiffness of arm-bridge angle (kT units)

asymmetry_stiffness = 100.0

Folding asymmetry stiffness of lower compartment (kT units)

average_steps_per_cycle() int
bridge_width = 7.5

Width of ATP bridge (nm)

cutoff4 = 3.5
cutoff5 = 3.5
cutoff6 = 3.0
cycles: int | None = 2

Number of SMC cycles (if set to None, will find approximate value using max_steps) Note: cycles are stochastic, so time per cycle is variable

diameter = 20

Diameter of initial loop (nm)

dna_config = 'folded'

configuration to generate

elbow_attraction = 30.0

Attractive energy between elbows in the APO state (kT units)

elbow_spacing = 2.5

Rest length between elbows in the APO state (nm)

elbows_stiffness = 30.0

Bending stiffness of elbows (kT units)

epsilon3 = 3.0
epsilon4 = 6.0
epsilon5 = 6.0
epsilon6 = 100.0
folding_angle_APO = 45.0

Folding angle of lower compartment (degrees)

folding_angle_ATP = 160.0

Folding angle of lower compartment (degrees)

folding_stiffness = 60.0

Folding stiffness of lower compartment (kT units)

force = 0.8

Stretching forces (pN) (set to any falsy value for no forces) WARNING: currently: if no forces -> ends are frozen

gamma = 0.5

Inverse of friction coefficient (ns)

hinge_radius = 1.5

Hinge radius (nm)

kB = 0.013806504

Boltzmann’s constant (pN nm / K)

kleisin_radius = 7.0

Radius of lower circular-arc compartment (nm)

loop = 100

Initial loop size (DNA beads)

max_steps: int | None = None

Max steps for run (None -> no maximum, will complete every cycle) Note: this is not a hard limit, some extra steps may be performed to complete a cycle

n = 5

Number of base pairs per DNA bead

non_random_steps = False

Disables the exponential sampling for APO,ATP,ADP steps

output_steps = 10000

Printing period (time steps)

rigid_hinge = True

hinge sections are connected by bonds

Type:

True

Type:

hinge is one rigid object, False

runs = 10

Number of independent runs

seed = 123

Random number seed

sigma = 2.5
sigma_SMC_DNA = 2.5

SMC-DNA hard-core repulsion radius = LJ sigma (nm)

site_cycle_period: int = 0

The number of SMC cycles between events where the cycling site is disabled. A value of zero disables this and uses the default site dynamics.

site_cycle_when: str = 'apo'

“apo”, “adp”.

Type:

When to re-enable the cycling site. Allowed values

site_stiffness = 100.0

Alignment stiffness of binding sites (kT units)

site_toggle_delay: int = 0

The number of SMC cycles between the cycling site being turned off and then on again. A value of zero means that the site will be enabled in the same cycle.

smc_force = 0.0

Extra force on SMC in the -x direction and +y direction (left & up)

spaced_beads_custom_stiffness: float = 1.0

multiple of the default DNA stiffness

spaced_beads_full_dna: bool = False

whether to place beads across the entire DNA length or not

spaced_beads_interval: int | None = None

how many DNA beads to leave between small obstacles

spaced_beads_size: float = 5.0

radius of beads along DNA (nm)

spaced_beads_smc_clearance: float = 5.0

length of bare DNA to keep next to SMC

spaced_beads_type = 1

rigid molecules

Type:

0

Type:

fene/expand bonds, 1

steps_ADP = 2000000

Average number of steps spent in ADP state (waiting for return to APO)

steps_APO = 2000000

Average number of steps spent in APO state (waiting for ATP binding)

steps_ATP = 8000000

Average number of steps spent in ATP state (waiting for ATP hydrolysis)

timestep = 0.0002

Simulation timestep (ns)

use_charges = False

Enable Coulomb interactions in LAMMPS

use_toroidal_hinge = True

old hinge type

Type:

True

Type:

toroidal hinge, False

smc_lammps.generate.generator module

generator.py

How to use:
  • create a Generator instance

  • create AtomType instances for each type you need

  • create BAI_Type for Bond/Angle/Improper types

  • define AtomGroup with positions and type

  • define BAI for interactions between specific atoms

  • define PairWise, and use add_interaction to define global interactions

  • append all defined AtomGroup, BAI, and PairWise instances to the lists in Generator

  • call Generator.write(file) to create the final datafile

class smc_lammps.generate.generator.AtomGroup(positions: ndarray[tuple[Any, ...], dtype[float32]], atom_type: AtomType, molecule_index: int, polymer_bond_type: BAI_Type | None = None, polymer_angle_type: BAI_Type | None = None, charge: float = 0.0)

Bases: object

Stores a list of atoms with the same AtomType.

Also supports polymer bond and angle for convenience.
positions

List of 3D atom positions (N, 3).

Type:

Nx3Array

atom_type

The atom type of all atoms in the group.

Type:

AtomType

molecule_index

The molecule index of all atoms in the group (see also Generator.molecule_override())

Type:

int

polymer_bond_type

BAI.Kind == BAI.BOND, forms a polymer in the order of the positions list.

Type:

BAI_Type | None

polymer_angle_type

BAI.Kind == BAI.ANGLE, adds angle potentials along the polymer.

Type:

BAI_Type | None

charge

The charge of all atoms in the group, used with atom style full.

Type:

float

__init__(positions: ndarray[tuple[Any, ...], dtype[float32]], atom_type: AtomType, molecule_index: int, polymer_bond_type: BAI_Type | None = None, polymer_angle_type: BAI_Type | None = None, charge: float = 0.0) None
property n: int

Number of atoms in the group.

Returns the length of the positions array.

Returns:

Number of atoms.

smc_lammps.generate.generator.AtomIdentifier

A unique identifier of an atom in a group.

alias of tuple[AtomGroup, int]

class smc_lammps.generate.generator.AtomType(mass: float = 1.0, unused: bool = False)

Bases: object

Represents a LAMMPS atom type.

See LAMMPS read_data ‘Atom Type Labels section’ and LAMMPS type labels.

exception UnusedIndex

Bases: AttributeError

LAMMPS does not allow gaps in atom types, therefore unused indices may cause issues.

__init__(mass: float = 1.0, unused: bool = False) None
property index: int

The LAMMPS atom type index.

If this is called for the first time, the global index is incremented.

Returns:

Index.

Raises:

self.UnusedIndex – Marked as unused, should not call this method.

class smc_lammps.generate.generator.BAI(type_: BAI_Type, *atoms: tuple[AtomGroup, int])

Bases: object

Represents a Bond/Angle/Improper interaction between a certain number of atoms.

See also PairWise for interactions between types of atoms (e.g. Lennard-Jones).

__init__(type_: BAI_Type, *atoms: tuple[AtomGroup, int]) None

Creates a Bond/Angle/Improper interaction between n atoms.

n depends on the BAI type:
  • Bond: n = 2

  • Angle: n = 3

  • Imprper: n = 4

Parameters:
  • type – The BAI type of the interaction.

  • atoms – List of n atoms.

atoms: list[tuple[AtomGroup, int]]
class smc_lammps.generate.generator.BAI_Kind(*values)

Bases: Enum

A LAMMPS Bond/Angle/Improper kind, see BAI_Type for usage.

ANGLE = 2

LAMMPS angle

BOND = 1

LAMMPS bond

IMPROPER = 3

LAMMPS improper (dihedral)

classmethod length_lookup(kind: BAI_Kind) int

Returns the number of arguments (atom ids) needed to define a BAI interaction.

Parameters:

kind – BAI Kind.

Returns:

Number of required arguments for the *_coeff command.

class smc_lammps.generate.generator.BAI_Type(kind: BAI_Kind, style: str, coefficients: str = '')

Bases: object

A LAMMPS Bond/Angle/Improper type.

_index

LAMMPS type index.

Type:

int

kind

Whether this is a bond, angle, or improper.

Type:

BAI_Kind

style

The style definition used in the LAMMPS command.

Type:

str

coefficients

The arguments to the LAMMPS command in style.

Type:

str

__init__(kind: BAI_Kind, style: str, coefficients: str = '') None

coefficients will be printed after the index in a datafile e.g. 3 harmonic 2 5 (where 3 is the index and the rest is the coefficients string)

get_string(omit_style: bool = False) str

Returns the string used with the *_coeff command.

Example:
>>> from smc_lammps.generate.generator import BAI_Type, BAI_Kind
>>> # create a harmonic bond with K=10.0 and a=1.0
>>> my_bond = BAI_Type(BAI_Kind.BOND, 'harmonic', '10.0 1.0')
>>> my_bond.get_string()
'1 harmonic 10.0 1.0'
Parameters:

omit_style – If True, do not print the BAI_Type.style.

Returns:

String in the ‘{index} {style} {coefficients}’ format.

property index: int

The LAMMPS Bond/Angle/Improper type index.

If this is called for the first time, the global index corresponding to the given BAI kind is incremented.

Returns:

Index.

indices = {BAI_Kind.ANGLE: 1, BAI_Kind.BOND: 1, BAI_Kind.IMPROPER: 1}

Global LAMMPS indices for each BAI kind.

smc_lammps.generate.generator.COORD_TYPE

The type to store coordinate positions.

class smc_lammps.generate.generator.Generator

Bases: object

Generates the LAMMPS data files from the given simulation properties.

class DynamicCoeffs(coeff_string: str, args: BAI_Type | list[AtomType])

Bases: object

Handles coefficients dynamically set within a LAMMPS script using pair_coeff, bond_coeff, angle_coeff, and improper_coeff commands.

__init__(coeff_string: str, args: BAI_Type | list[AtomType]) None
classmethod create_from_pairwise(pairwise: PairWise, type1: AtomType, type2: AtomType, values: list[Any] | None) DynamicCoeffs
static get_script_bai_command_name(pair_or_BAI: BAI_Kind | None) str

Returns the command name for a pair / BAI interaction. Used to redefine coefficients within a LAMMPS script.

write_script_bai_coeffs(file: TextIO) None

Writes a LAMMPS command for a pair / BAI interaction to a file.

__init__() None
add_atom_groups(*args: AtomGroup) None

Adds new atom groups.

Atom groups with empty positions are ignored.

atom_group_map: list[int]

Maps the index of an atom group to the global index offset for the LAMMPS atom indices. Note! This is only defined after calling Generator.write_atoms().

bais: list[BAI]

List of Bond/Angle/Improper interactions in the simulation.

box_width

Size of the simulation box (cube side length).

charge_override: dict[tuple[AtomGroup, int], float]

Individual charge value overrides for atoms. Takes precedence over the AtomGroup charge value.

check_charges() None

Checks for discrepancies between the use_charges flag and the actual atom charges.

get_BAI_coeffs_header(kind: BAI_Kind) str

Returns the header string corresponding to a Bond/Angle/Improper kind.

static get_BAI_header(kind: BAI_Kind) str
static get_BAI_style_command_name(kind: BAI_Kind) str
get_BAI_styles_command() str

Returns the command to define the BAI styles in a LAMMPS script.

Returns:

LAMMPS command string.

get_all_BAI_styles() dict[BAI_Kind, list[str]]

Returns a list of unique styles for each BAI kind.

Returns:

Dictionary which maps each BAI kind to a list of unique styles.

get_all_atom_types() list[AtomType]

Returns a list of all atom types across all atom groups.

Returns:

List of atom types, sorted by AtomType.index.

get_all_types(kind: BAI_Kind) list[BAI_Type]

Returns a list of all Bond/Angle/Improper types across all atom groups.

Parameters:

kind – BAI kind to filter by.

Returns:

List of BAI types of the given kind, sorted by BAI_Type.index.

get_amounts() tuple[int, int, int]

Returns the total amount of bonds, angles, and impropers.

Returns:

Tuple of (# bonds, # angles, # impropers)

get_atom_index(atom_id: tuple[AtomGroup, int]) int

Returns the absolute LAMMPS index for an atom.

Attention

You must call write_atoms() before using this method.

get_atom_style() str

Returns the atom_style for LAMMPS (based on use_charges).

get_atom_style_command() str
get_atoms_header() str
get_bai_dict_by_type() dict[BAI_Kind, list[BAI]]

Returns a dictionary mapping the Bond/Angle/Improper kind to a list of all BAIs which have that type.

Returns:

Dictionary which maps kind to BAIs of that kind.

get_hybrid_or_single_style() dict[BAI_Kind, str]

Returns the style needed for the BAI Coeffs header.

get_total_atoms() int

Returns the total number of atoms across all groups.

hybrid_styles

Used when multiple styles are defined, default of ‘hybrid’ should work in most cases.

molecule_override: dict[tuple[AtomGroup, int], int]

Individual molecule id overrides for atoms. Takes precedence over the AtomGroup molecule id.

move_all_atoms(shift: ndarray[tuple[Any, ...], dtype[float32]]) None

Moves all atoms in a certain direction.

Useful to place the system at the center of the simulation box.

Parameters:

shift – All atoms are shifted by this vector (1, 3).

pair_interactions: list[PairWise]

List of pair interactions in the simulation.

random_shift

Function that returns a random shift vector. This is useful to avoid exact overlap, which causes LAMMPS to crash during kspace calculations (e.g. with pair_style coul).

set_system_size(box_width: float) None

Sets the box size of the simulation.

Note: box is cubic.

Parameters:

box_width – Width of box.

Raises:

ValueError – Received negative or zero box_width.

use_charges

Whether to enable charges (atom_style ‘full’) or not (atom_style ‘molecule’).

write_BAI_coeffs(file: TextIO) None

Writes the Bond/Angle/Improper coefficients for each BAI kind to a file.

write_amounts(file: TextIO) None

Writes the amount of atoms, bonds, angles, and impropers to a file.

Parameters:

file – File to write to.

write_atoms(file: TextIO) None

Writes the Atoms header and all atom positions to a file.

Attention

This calls _set_up_atom_group_map(), read the note there.

Parameters:

file – File to write to.

write_bai(file: TextIO) None

Writes the Bond/Angle/Improper headers and all corresponding BAI interactions to a file.

write_coeffs(file: TextIO) None

Writes the coefficient information to a file.

Useful when restarting a simulation.

write_full(file: TextIO) None

Writes a full LAMMPS data file.

write_header(file: TextIO) None

Writes the top header to a file.

Note: LAMMPS always ignores the first line of a data file, which should be this header.

Parameters:

file – File to write to.

write_masses(file: TextIO) None

Writes the masses of all atom types to a file.

Parameters:

file – File to write to.

write_pair_interactions(file: TextIO) None

Writes the Pair Coeffs header(s) and corresponding pair interactions to a file.

write_positions_and_bonds(file: TextIO) None

Writes the positions and bonds to a file.

write_system_size(file: TextIO) None

Writes the system size to a file.

Parameters:

file – File to write to.

Raises:

TypeError – The box_width is None.

write_types(file: TextIO) None

Writes the amount of atom types, bond types, angle types, and improper types to a file.

Parameters:

file – File to write to.

class smc_lammps.generate.generator.MoleculeId

Bases: object

Represents a LAMMPS molecule id.

See LAMMPS molecule.

classmethod get_next() int

Increments the global index and returns it.

Returns:

Index.

smc_lammps.generate.generator.Nx3Array

An (N, 3) array of positions.

alias of ndarray[tuple[Any, …], dtype[float32]]

class smc_lammps.generate.generator.PairWise(header: str, template: str, default: list[Any] | None)

Bases: object

Represents pair interactions between all atoms of two atoms ids.

header

Definition of the interaction style, e.g. 'PairIJ Coeffs # hybrid'.

Type:

str

template

Format string with empty formatters {} for the interaction parameters, e.g. 'lj/cut {} {} {}'.

Type:

str

default

List of default parameters. If None, do not insert missing interactions. This is used to fill out all interactions, since LAMMPS requires them to all be explicitly defined.

Type:

list[Any] | None

__init__(header: str, template: str, default: list[Any] | None) None
add_interaction(atom_type1: AtomType, atom_type2: AtomType, *args: Any, **kwargs: bool) PairWise

Adds an iteraction.

Indices are sorted automatically, which is required by LAMMPS.

get_all_interaction_pairs(all_atom_types: list[AtomType]) list[tuple[AtomType, AtomType]]

Returns all possible interactions, whether they are set by the user or not.

Parameters:

all_atom_types – List of all the atom types that exist in the simulation.

Returns:

List of interactions (t_1, t_2) with t_1 <= t_2.

get_all_interactions(all_atom_types: list[AtomType]) list[tuple[AtomType, AtomType, str]]

Returns actual interactions to define.

Applies the default where no interaction was specified by the user.

Parameters:

all_atom_types – List of all the atom types that exist in the simulation.

Returns:

List of (t_1, t_2, f) where t_1 <= t_2 and f is the formatted string.

pair_in_inter(interaction: tuple[AtomType, AtomType]) tuple[AtomType, AtomType, list[Any]] | None

Checks if an interaction is defined in pairs.

Parameters:

interaction – Interaction to look for.

Returns:

The pair as it is defined in pairs (may have opposite order), or None if no pair was found.

pairs: list[tuple[AtomType, AtomType, list[Any]]]
write(file: TextIO, atom_types: list[AtomType]) None

Writes the Pair Coeffs header and all pair interactions to a file.

Example:

PairIJ Coeffs 1 8 lj/cut 0.0 0.0 0.0

smc_lammps.generate.generator.all_tests()
smc_lammps.generate.generator.main()
smc_lammps.generate.generator.test_simple_atoms()
smc_lammps.generate.generator.test_simple_atoms_polymer()
smc_lammps.generate.generator.test_with_bonds()
smc_lammps.generate.generator.test_with_pairs()
smc_lammps.generate.generator.write_if_non_zero(file: TextIO, fmt_string: str, amount: int)

Writes an amount to a file only if it is not zero.

Parameters:
  • file – File to write to.

  • fmt_string – Format string with one empty formatter {} for the amount.

  • amount – The amount to write if non-zero.

smc_lammps.generate.generator.Generator._set_up_atom_group_map(self) None

Sets the atom_group_map based on the current atom groups.

Attention

The atom groups must not change after calling this function.

smc_lammps.generate.parameters_template module

smc_lammps.generate.util module

smc_lammps.generate.util.create_phase(phase_path: Path, options: Sequence[DynamicCoeffs])

creates a file containing coefficients to dynamically load in LAMMPS scripts

smc_lammps.generate.util.create_phase_wrapper(phase_path: Path, options: Sequence[DynamicCoeffs | None])

filters out None values and then calls create_phase

smc_lammps.generate.util.get_closest(array, position) int

returns the index of the array that is closest to the given position

smc_lammps.generate.util.get_parameters(path: Path) Parameters

Load the parameters from a parameters.py file.

smc_lammps.generate.util.get_project_root() Path
smc_lammps.generate.util.pos_from_id(atom_id: tuple[AtomGroup, int]) ndarray[tuple[Any, ...], dtype[float32]]

get the position of an atom from its identifier