reforge.forge package

Submodules

reforge.forge.cgmap module

CG mapping tools

Description:

Provides functions to map atomistic PDB data to a coarse-grained (CG) representation.

Author: DY

reforge.forge.cgmap.move_o3(system)[source]

Move each O3’ atom to the next residue. Needed for some CG nucleic force fields due to the phosphate group mapping.

Parameters:

system: A system object containing chains and residues.

reforge.forge.cgmap.map_residue(residue, mapping, atid)[source]

Map an atomistic residue to a coarse-grained (CG) residue.

For each bead defined in the mapping dictionary, this function creates a new bead atom. The bead’s coordinates are determined by averaging the coordinates of the atoms in the original residue that match the bead’s atom names.

Parameters:

residue: The original residue object. mapping (dict): Dictionary with bead names as keys and lists of atom names as values. atid (int): Starting atom id to assign to new beads.

Returns:

list: List of new bead atoms representing the coarse-grained residue.

reforge.forge.cgmap.map_chain(chain, ff, atid=1)[source]

Map a chain of atomistic residues to a coarse-grained (CG) representation.

For each residue in the chain, the function retrieves the corresponding bead mapping from the force field (ff) based on the residue’s name. For the first residue, the mapping for “BB1” is removed. Each residue is then converted to its CG representation using the map_residue function, and the resulting beads are collected into a single list.

Parameters:

chain: List of residue objects. ff: Force field object that contains a ‘mapping’ dictionary keyed by residue name. atid (int): Starting atom id for new beads (default is 1).

Returns:

list: List of coarse-grained bead atoms representing the entire chain.

reforge.forge.cgmap.map_model(model, ff, atid=1)[source]

Map a model of atomistic residues to a coarse-grained (CG) representation.

Parameters:

model: Iterable of chains (each a list of residue objects). ff: Force field object containing a ‘mapping’ dictionary keyed by residue name. atid (int): Starting atom id for new beads (default is 1).

Returns:

list: List of coarse-grained bead atoms representing the entire model.

reforge.forge.forcefields module

Collection of CG force fields

Description:

This module defines force fields for Martini RNA and nucleic acids. It provides classes for reading ITP files and organizing force-field parameters for coarse-grained simulations.

Author: DY

reforge.forge.forcefields.nsplit(*x)[source]
class reforge.forge.forcefields.NucleicForceField(directory, mol, version)[source]

Bases: object

Base class for nucleic acid force fields.

static read_itp(resname, directory, mol, version)[source]

Read an ITP file for a given residue.

Args:

resname (str): Residue name. directory (str): Directory where force field files are stored. mol (str): Molecule name. version (str): Version string.

Returns:

dict: Parsed ITP data.

static read_itps(mol, ff_type, version)[source]

Read ITP files for nucleobases A, C, G, and U.

Args:

mol (str): Molecule name. ff_type (str): Force-field type (e.g. “polar”). version (str): Version string.

Returns:

tuple: (a_itp, c_itp, g_itp, u_itp)

static itp_to_indata(itp_data)[source]

Convert ITP data into individual parameter lists.

Args:

itp_data (dict): Parsed ITP file data.

Returns:

tuple: (sc_bonds, sc_angles, sc_dihs, sc_cons, sc_excls, sc_pairs, sc_vs3s)

static parameters_by_resname(resnames, directory, mol, version)[source]

Obtain parameters for each residue.

Args:

resnames (list): List of residue names. directory (str): Directory containing ITP files. mol (str): Molecule name. version (str): Version string.

Returns:

dict: Mapping of residue name to its parameters.

sc_bonds(resname)[source]
sc_angles(resname)[source]
sc_dihs(resname)[source]
sc_cons(resname)[source]
sc_excls(resname)[source]
sc_pairs(resname)[source]
sc_vs3s(resname)[source]
sc_blist(resname)[source]

Get all bonded parameters for a residue.

Args:

resname (str): Residue name.

Returns:

list: List of bonded parameter lists.

class reforge.forge.forcefields.Martini30RNA(directory='rna_reg', mol='test', version='new')[source]

Bases: NucleicForceField

Force field for Martini 3.0 RNA.

resnames = ['A', 'C', 'G', 'U']
alt_resnames = ['ADE', 'CYT', 'GUA', 'URA']
bb_mapping = {'BB1': ('P', 'OP1', 'OP2', "O5'", "O3'", 'O1P', 'O2P'), 'BB2': ("C5'", "1H5'", "2H5'", "H5'", "H5''", "C4'", "H4'", "O4'", "C3'", "H3'"), 'BB3': ("C1'", "C2'", "O2'", "O4'")}
a_mapping = {'SC1': ('N9', 'C8', 'H8'), 'SC2': ('N3', 'C4'), 'SC3': ('N1', 'C2', 'H2'), 'SC4': ('N6', 'C6', 'H61', 'H62'), 'SC5': ('N7', 'C5')}
c_mapping = {'SC1': ('N1', 'C5', 'C6'), 'SC2': ('C2', 'O2'), 'SC3': ('N3',), 'SC4': ('N4', 'C4', 'H41', 'H42')}
g_mapping = {'SC1': ('C8', 'H8', 'N9'), 'SC2': ('C4', 'N3'), 'SC3': ('C2', 'N2', 'H21', 'H22'), 'SC4': ('N1',), 'SC5': ('C6', 'O6'), 'SC6': ('C5', 'N7')}
u_mapping = {'SC1': ('N1', 'C5', 'C6'), 'SC2': ('C2', 'O2'), 'SC3': ('N3',), 'SC4': ('C4', 'O4')}
mapping = {'A': {'BB1': ('P', 'OP1', 'OP2', "O5'", "O3'", 'O1P', 'O2P'), 'BB2': ("C5'", "1H5'", "2H5'", "H5'", "H5''", "C4'", "H4'", "O4'", "C3'", "H3'"), 'BB3': ("C1'", "C2'", "O2'", "O4'"), 'SC1': ('N9', 'C8', 'H8'), 'SC2': ('N3', 'C4'), 'SC3': ('N1', 'C2', 'H2'), 'SC4': ('N6', 'C6', 'H61', 'H62'), 'SC5': ('N7', 'C5')}, 'ADE': {'BB1': ('P', 'OP1', 'OP2', "O5'", "O3'", 'O1P', 'O2P'), 'BB2': ("C5'", "1H5'", "2H5'", "H5'", "H5''", "C4'", "H4'", "O4'", "C3'", "H3'"), 'BB3': ("C1'", "C2'", "O2'", "O4'"), 'SC1': ('N9', 'C8', 'H8'), 'SC2': ('N3', 'C4'), 'SC3': ('N1', 'C2', 'H2'), 'SC4': ('N6', 'C6', 'H61', 'H62'), 'SC5': ('N7', 'C5')}, 'C': {'BB1': ('P', 'OP1', 'OP2', "O5'", "O3'", 'O1P', 'O2P'), 'BB2': ("C5'", "1H5'", "2H5'", "H5'", "H5''", "C4'", "H4'", "O4'", "C3'", "H3'"), 'BB3': ("C1'", "C2'", "O2'", "O4'"), 'SC1': ('N1', 'C5', 'C6'), 'SC2': ('C2', 'O2'), 'SC3': ('N3',), 'SC4': ('N4', 'C4', 'H41', 'H42')}, 'CYT': {'BB1': ('P', 'OP1', 'OP2', "O5'", "O3'", 'O1P', 'O2P'), 'BB2': ("C5'", "1H5'", "2H5'", "H5'", "H5''", "C4'", "H4'", "O4'", "C3'", "H3'"), 'BB3': ("C1'", "C2'", "O2'", "O4'"), 'SC1': ('N1', 'C5', 'C6'), 'SC2': ('C2', 'O2'), 'SC3': ('N3',), 'SC4': ('N4', 'C4', 'H41', 'H42')}, 'G': {'BB1': ('P', 'OP1', 'OP2', "O5'", "O3'", 'O1P', 'O2P'), 'BB2': ("C5'", "1H5'", "2H5'", "H5'", "H5''", "C4'", "H4'", "O4'", "C3'", "H3'"), 'BB3': ("C1'", "C2'", "O2'", "O4'"), 'SC1': ('C8', 'H8', 'N9'), 'SC2': ('C4', 'N3'), 'SC3': ('C2', 'N2', 'H21', 'H22'), 'SC4': ('N1',), 'SC5': ('C6', 'O6'), 'SC6': ('C5', 'N7')}, 'GUA': {'BB1': ('P', 'OP1', 'OP2', "O5'", "O3'", 'O1P', 'O2P'), 'BB2': ("C5'", "1H5'", "2H5'", "H5'", "H5''", "C4'", "H4'", "O4'", "C3'", "H3'"), 'BB3': ("C1'", "C2'", "O2'", "O4'"), 'SC1': ('C8', 'H8', 'N9'), 'SC2': ('C4', 'N3'), 'SC3': ('C2', 'N2', 'H21', 'H22'), 'SC4': ('N1',), 'SC5': ('C6', 'O6'), 'SC6': ('C5', 'N7')}, 'U': {'BB1': ('P', 'OP1', 'OP2', "O5'", "O3'", 'O1P', 'O2P'), 'BB2': ("C5'", "1H5'", "2H5'", "H5'", "H5''", "C4'", "H4'", "O4'", "C3'", "H3'"), 'BB3': ("C1'", "C2'", "O2'", "O4'"), 'SC1': ('N1', 'C5', 'C6'), 'SC2': ('C2', 'O2'), 'SC3': ('N3',), 'SC4': ('C4', 'O4')}, 'URA': {'BB1': ('P', 'OP1', 'OP2', "O5'", "O3'", 'O1P', 'O2P'), 'BB2': ("C5'", "1H5'", "2H5'", "H5'", "H5''", "C4'", "H4'", "O4'", "C3'", "H3'"), 'BB3': ("C1'", "C2'", "O2'", "O4'"), 'SC1': ('N1', 'C5', 'C6'), 'SC2': ('C2', 'O2'), 'SC3': ('N3',), 'SC4': ('C4', 'O4')}}
sc_atoms(resname)[source]

Return side-chain atoms for the given residue.

class reforge.forge.forcefields.Martini31Nucleic[source]

Bases: NucleicForceField

Force field for Martini 3.1 nucleic acids.

bb_mapping = [['P', 'OP1', 'OP2', "O5'", "O3'", 'O1P', 'O2P'], ["C5'", "1H5'", "2H5'", "H5'", "H5''", "C4'", "H4'", "O4'", "C3'", "H3'"], ["C1'", "C2'", "O2'", "O4'"]]
mapping = {'A': {'SC1': ['TA1', 'TA2', 'TA3', 'TA4', 'TA5', 'TA6']}, 'C': {'SC1': ['TY1', 'TY2', 'TY3', 'TY4', 'TY5']}, 'G': {'SC1': ['TG1', 'TG2', 'TG3', 'TG4', 'TG5', 'TG6', 'TG7', 'TG8']}, 'U': {'SC1': ['TU1', 'TU2', 'TU3', 'TU4', 'TU5', 'TU6', 'TU7']}}
sc_atoms(resname)[source]

Return side-chain atoms for a given residue.

reforge.forge.geometry module

Description:

Provides geometric calculations (distances, angles, dihedrals) for coarse-grained systems and functions to compute bond lists from a system (using a reference topology).

Author: DY

reforge.forge.geometry.get_distance(v1, v2)[source]

Compute the Euclidean distance between two vectors (returned in nanometers).

Args:

v1: First vector (list or array-like). v2: Second vector (list or array-like).

Returns:

float: Distance between v1 and v2 divided by 10.

reforge.forge.geometry.get_angle(v1, v2, v3)[source]

Calculate the angle at vertex v2 (in degrees) given three points.

Args:

v1: First vector (list or array-like). v2: Vertex vector (list or array-like). v3: Third vector (list or array-like).

Returns:

float: Angle in degrees.

reforge.forge.geometry.get_dihedral(v1, v2, v3, v4)[source]

Calculate the dihedral angle (in degrees) defined by four points.

Args:

v1: First vector. v2: Second vector. v3: Third vector. v4: Fourth vector.

Returns:

float: Dihedral angle in degrees.

reforge.forge.geometry.calc_bonds(atoms, bonds)[source]

Calculate bond distances from a topology bonds object and a list of atoms.

Args:

atoms: List of atom objects. bonds: Topology bonds object with attributes ‘conns’, ‘params’, and ‘comms’.

Returns:

BondList: A list of bonds with computed distances.

reforge.forge.geometry.calc_angles(atoms, angles)[source]

Calculate bond angles from a topology angles object and a list of atoms.

Args:

atoms: List of atom objects. angles: Topology angles object with attributes ‘conns’, ‘params’, and ‘comms’.

Returns:

BondList: A list of angles with computed values.

reforge.forge.geometry.calc_dihedrals(atoms, dihs)[source]

Calculate dihedral angles from a topology dihedrals object and a list of atoms.

Args:

atoms: List of atom objects. dihs: Topology dihedrals object with attributes ‘conns’, ‘params’, and ‘comms’.

Returns:

BondList: A list of dihedrals with computed values.

reforge.forge.geometry.get_cg_bonds(inpdb, top)[source]

Calculate bonds, angles, and dihedrals for a coarse-grained system.

Args:

inpdb (str): Input PDB file for the coarse-grained system. top: Topology object containing bond, angle, and dihedral definitions.

Returns:

tuple: Three BondList objects: (bonds, angles, dihedrals).

reforge.forge.geometry.get_aa_bonds(inpdb, ff, top)[source]

Calculate bonds, angles, and dihedrals for an all-atom system.

Args:

inpdb (str): Input PDB file for the all-atom system. ff: Force field object. top: Topology object containing reference bonds, angles, and dihedrals.

Returns:

tuple: Three BondList objects: (bonds, angles, dihedrals).

reforge.forge.persistence_length module

Description:

This module computes the persistence length of a nucleic acid system by calculating angles between consecutive bond vectors in a coarse-grained model. It provides functions to read a structure from a PDB file, compute geometric measures (angles, rotation matrices), and finally plot the average angles.

Author: DY

reforge.forge.persistence_length.get_structure_cif(cif_id='path/to/cif')[source]

Read a structure from a CIF file.

Args:

cif_id (str): Path to the CIF file.

Returns:

structure: A Bio.PDB structure object.

reforge.forge.persistence_length.get_structure_pdb(pdb_id='path/to/pdb')[source]

Read a structure from a PDB file.

Args:

pdb_id (str): Path to the PDB file.

Returns:

structure: A Bio.PDB structure object.

reforge.forge.persistence_length.get_residues(model)[source]

Extract all residues from a model.

Args:

model: A Bio.PDB model object.

Returns:

list: A list of residue objects.

reforge.forge.persistence_length.get_resid(residue)[source]

Extract the residue id from the string representation of a residue.

Args:

residue: A Bio.PDB residue object.

Returns:

int: The residue id.

reforge.forge.persistence_length.get_atoms_by_name(residues, atom_name='BB3')[source]

Select atoms with a given name from a list of residues.

Args:

residues (list): List of residue objects. atom_name (str): Name of the atom to select (default “BB3”).

Returns:
tuple: A tuple (atoms, resids) where atoms is a list of matching atom objects,

and resids is a list of corresponding residue ids.

reforge.forge.persistence_length.get_coords(atoms)[source]

Get the coordinates for a list of atoms.

Args:

atoms (list): List of atom objects.

Returns:

list: List of coordinate arrays.

reforge.forge.persistence_length.get_angle(v1, v2)[source]

Compute the angle (in degrees) between two vectors.

Args:

v1: First vector. v2: Second vector.

Returns:

float: Angle in degrees.

reforge.forge.persistence_length.rotation_matrix_from_vectors(vec1, vec2)[source]

Find the rotation matrix that aligns vec1 to vec2.

Args:

vec1 (array-like): Source vector. vec2 (array-like): Destination vector.

Returns:

np.ndarray: A 3x3 rotation matrix.

reforge.forge.persistence_length.persistence_length(structure)[source]

Calculate and plot the persistence length from a structure.

For each model in the structure, the function computes the angle between bond vectors (taken from atoms with name “BB3”) after aligning them with a reference vector. The average angles are then plotted and saved to a file.

Args:

structure: A Bio.PDB structure object.

reforge.forge.persistence_length.test_rmats()[source]

Test function for rotation_matrix_from_vectors.

reforge.forge.persistence_length.main()[source]

Main function to compute persistence length.

reforge.forge.topology module

Topology Module

Description:

This module provides classes and functions to construct a coarse-grained topology from force field parameters. It defines the Topology class along with a helper BondList class. Methods are provided to process atoms, bonds, and connectivity information, and to generate topology files for coarse-grained simulations.

Usage Example:
>>> from topology import Topology
>>> from reforge.forcefield import NucleicForceField
>>> ff = NucleicForceField(...)  # Initialize the force field instance
>>> topo = Topology(ff, sequence=['A', 'T', 'G', 'C'])
>>> topo.from_sequence(['A', 'T', 'G', 'C'])
>>> topo.write_to_itp('output.itp')
Requirements:
  • Python 3.x

  • NumPy

  • reForge utilities and force field modules

Author: DY Date: YYYY-MM-DD

class reforge.forge.topology.BondList(iterable=(), /)[source]

Bases: list

BondList Class

Description:

A helper subclass of the built-in list for storing bond information. Each element represents a bond in the form [connectivity, parameters, comment].

Variables:
  • list) ((inherited from) – Holds individual bond representations.

  • Example (Usage) –

    >>> bonds = BondList([['C1-C2', [1.0, 1.5], 'res1 bead1'],
    ...                   ['C2-O1', [1.1, 1.6], 'res2 bead2']])
    >>> print(bonds.conns)
    ['C1-C2', 'C2-O1']
    >>> bonds.conns = ['C1-C2_mod', 'C2-O1_mod']
    >>> print(bonds.conns)
    ['C1-C2_mod', 'C2-O1_mod']
    

property conns

Get connectivity values from each bond.

Returns:

List of connectivity values (index 0 of each bond).

Return type:

list

property params

Get parameters from each bond.

Returns:

List of parameter values (index 1 of each bond).

Return type:

list

property comms

Get comments from each bond.

Returns:

List of comments (index 2 of each bond).

Return type:

list

property measures

Get measure values from each bond.

Returns:

List of measures extracted from the second element of the bond parameters.

Return type:

list

categorize()[source]

Categorize bonds based on their comments.

Description:

Uses the stripped comment (index 2) as a key to group bonds into a dictionary.

Returns:

A dictionary mapping each unique comment to a BondList of bonds.

Return type:

dict

filter(condition, bycomm=True)[source]

Filter bonds based on a provided condition.

Description:

Selects bonds for which the condition (a callable) returns True. By default, the condition is applied to the comment field.

Parameters:
  • condition (callable) – A function that takes a bond (or its comment) as input and returns True if the bond should be included.

  • bycomm (bool, optional) – If True, the condition is applied to the comment (default is True).

Returns:

A new BondList containing the bonds that meet the condition.

Return type:

BondList

class reforge.forge.topology.Topology(forcefield, sequence: List = None, secstruct: List = None, **kwargs)[source]

Bases: object

Topology Class

Description:

Constructs a coarse-grained topology from force field parameters. Provides methods for processing atoms, bonds, and connectivity, and for generating a topology file for coarse-grained simulations.

Variables:
  • ff (object) – Force field instance.

  • sequence (list) – List of residue names.

  • name (str) – Molecule name.

  • nrexcl (int) – Exclusion parameter.

  • atoms (list) – List of atom records.

  • elnet (bonds, angles, dihs, cons, excls, pairs, vs3s, posres,) – BondList instances for various bonded interactions.

  • blist (list) – List containing all bond-type BondLists.

  • secstruct (list) – Secondary structure as a list of characters.

  • natoms (int) – Total number of atoms.

lines() list[source]

Generate the topology file as a list of lines.

Returns:

A list of strings, each representing a line in the topology file.

Return type:

list

write_to_itp(filename: str)[source]

Write the topology to an ITP file.

Parameters:

filename (str) – The output file path.

process_atoms(start_atom: int = 0, start_resid: int = 1)[source]

Process atoms based on the sequence and force field definitions.

Description:

For each residue in the sequence, constructs atom records using both backbone and sidechain definitions from the force field.

Parameters:
  • start_atom (int, optional) – Starting atom ID (default is 0).

  • start_resid (int, optional) – Starting residue ID (default is 1).

process_bb_bonds(start_atom: int = 0, start_resid: int = 1)[source]

Process backbone bonds using force field definitions.

Parameters:
  • start_atom (int, optional) – Starting atom ID.

  • start_resid (int, optional) – Starting residue ID.

process_sc_bonds(start_atom: int = 0, start_resid: int = 1)[source]

Process sidechain bonds using force field definitions.

Parameters:
  • start_atom (int, optional) – Starting atom ID.

  • start_resid (int, optional) – Starting residue ID.

elastic_network(atoms, anames: List[str] = None, el: float = 0.5, eu: float = 1.1, ef: float = 500)[source]

Construct an elastic network between selected atoms.

Parameters:
  • atoms (list) – List of PDB atom objects.

  • anames (list of str, optional) – Atom names to include (default: [“BB1”, “BB3”]).

  • el (float, optional) – Lower distance cutoff.

  • eu (float, optional) – Upper distance cutoff.

  • ef (float, optional) – Force constant.

from_sequence(sequence, secstruc=None)[source]

Build topology from a given sequence.

Parameters:
  • sequence (list) – Nucleic acid sequence (list of residue names).

  • secstruc (list, optional) – Secondary structure. If not provided, defaults to all ‘F’.

from_chain(chain, secstruc=None)[source]

Build topology from a chain instance.

Parameters:
  • chain (object) – Chain object containing residues.

  • secstruc (list, optional) – Secondary structure.

static merge_topologies(topologies)[source]

Merge multiple Topology instances into one.

Parameters:

topologies (list) – List of Topology objects.

Returns:

A merged Topology instance.

Return type:

Topology

Module contents