mchem

Core

CLI entry point: convert PDB to SQLite DB, solvate PDB, etc.

System of force terms and particles; load/save from SQLite DB.

class mchem.system.Box(a: float, b: float, c: float, alpha: float, beta: float, gamma: float)[source]

Bases: object

Periodic box dimensions from PDB CRYST1: cell lengths a, b, c (Angstroms) and angles alpha, beta, gamma (degrees).

a: float
b: float
c: float
alpha: float
beta: float
gamma: float
__init__(a: float, b: float, c: float, alpha: float, beta: float, gamma: float) None
mchem.system.register_term_class(cls)[source]

Register a custom term dataclass so it can be deserialized from a database file. Call this before loading a .db that contains the class.

class mchem.system.System(path: PathLike | None = None)[source]

Bases: object

Container for force-field terms and metadata, backed by SQLite.

Terms are stored by table name (e.g. Particle, AmoebaBond). Metadata is stored in a separate meta table. Use register_term_class() to register custom term dataclasses before loading a DB that contains them.

__init__(path: PathLike | None = None)[source]

Create an empty system or load from an existing SQLite DB.

Parameters:

path (os.PathLike, optional) – If given, path to SQLite file to load; otherwise an empty system.

__getitem__(key: str)[source]

Return the term list for the given table name.

property data

Mapping of table name to TermList of terms.

property meta

Mapping of metadata keys to values (e.g. name, units).

fromDatabase()[source]

Load all tables from the connected SQLite database into data and meta.

getTableNames()[source]

Return list of table names in the connected database.

getTermsByName(name: str) TermList[source]

Return the TermList for the given term table name.

addTerm(term, name: str | None = None)[source]

Append a single term; use its class name as table name if name is not given.

addTerms(terms: TermList, name: str | None = None)[source]

Append all terms from a TermList.

addMeta(key: str, value: Any)[source]

Set a metadata key (hyphens in key are replaced with underscores).

save(path: PathLike, overwrite: bool = False)[source]

Write meta and all term tables to a new SQLite file.

Parameters:
  • path (os.PathLike) – Output DB path.

  • overwrite (bool, optional) – If True, replace existing file; otherwise raise if file exists.

close()[source]

Close the database connection.

Residue templates (atoms and bonds) loaded from XML; used for topology matching.

class mchem.template.ResidueTemplate(name: str, atoms: Dict[str, Dict[str, Any]], bonds: List[Dict[str, Any]], altNames: List[str] = [])[source]

Bases: object

Template for a residue: atom names (with alternates) and bond list for matching.

__init__(name: str, atoms: Dict[str, Dict[str, Any]], bonds: List[Dict[str, Any]], altNames: List[str] = [])[source]
Parameters:
  • name (str) – Residue template name.

  • atoms (dict) – Map atom name -> dict with altNames (and optionally atomType).

  • bonds (list) – List of dicts with keys atom1, atom2, order.

  • altNames (list, optional) – Alternative names for this template.

classmethod fromXMLData(data: Element)[source]

Build a ResidueTemplate from an XML Residue element (name, Atom, Bond children).

setAtomType(atomName: str, atomType: str)[source]

Set force-field atom type for the given atom name.

getAtomType(atomName: str)[source]

Return force-field atom type for the given atom name.

mchem.template.loadTemplateDefinitions(fname: PathLike)[source]

Load residue templates from an XML file and add them to TEMPLATES.

Parameters:

fname (os.PathLike) – Path to XML file containing Residue elements.

Molecular topology: residues, atoms, bonds, and connectivity.

class mchem.topology.Residue(name: str, number: int, chain: str = 'A', insertionCode: str = '')[source]

Bases: object

A residue (e.g. amino acid) containing atoms and optional link to a Topology.

__init__(name: str, number: int, chain: str = 'A', insertionCode: str = '')[source]
Parameters:
  • name (str) – Residue name (e.g. "ALA").

  • number (int) – Residue sequence number.

  • chain (str, optional) – Chain identifier.

  • insertionCode (str, optional) – PDB insertion code.

property name: str
setName(newName: str)[source]
property stdName: str
setStdName(newName: str)[source]
property topology
setTopology(topology)[source]
hasTopology()[source]
property idx: int
setIdx(newIdx: int)[source]
addAtom(atom)[source]
getAtomWithName(name: str)[source]
removeAtom(atom)[source]
property numAtoms
property atoms
matchTemplate(template: ResidueTemplate, stdResidueName: bool = True)[source]

Match this residue to a template: align atom names (including alternates) and add template bonds to topology.

Parameters:
  • template (ResidueTemplate) – Template to match (atom names and bonds).

  • stdResidueName (bool, optional) – If True, set residue name to template name.

Returns:

True if match succeeded.

Return type:

bool

class mchem.topology.Atom(name: str, element: str)[source]

Bases: object

An atom with element, position, optional residue/topology, bonds, and force-field type/class.

__init__(name: str, element: str)[source]
Parameters:
  • name (str) – Atom name (e.g. "CA").

  • element (str) – Element symbol (key in ELEMENTS).

setPosition(pos)[source]
property polarizationGroup
setPolarizationGroup(setOfAtoms)[source]
addToPolarizationGroup(other)[source]
property mass: float
property symbol: str
property atomicNum: int
property atomType
property atomClass
setAtomType(atype)[source]
setAtomClass(aclass)[source]
property parametrized
property name: str
setName(newname: str)[source]
property idx: int
setIdx(newIdx: int)[source]
property residue
hasResidue()[source]
setResidue(residue: None | Residue)[source]
property topology
hasTopology()[source]
property bonds
generateNeighbors()[source]
getNeighbors()[source]
getHighOrderNeighbors(order: int)[source]
generateTopologicalInfo(maxConnect: int)[source]
property pathsToBondedAtoms
property bondedAtoms
getTopologicalInfo()[source]
addBond(bond)[source]
pruneBonds()[source]
class mchem.topology.Bond(atom1: Atom, atom2: Atom, order: float)[source]

Bases: object

A bond between two atoms with a bond order.

__init__(atom1: Atom, atom2: Atom, order: float)[source]
Parameters:
  • atom1 (Atom) – Bonded atoms.

  • atom2 (Atom) – Bonded atoms.

  • order (float) – Bond order (e.g. 1.0, 2.0).

hasTopology() bool[source]
mchem.topology.require_editable(func)[source]

Decorator: raise if topology is not editable.

mchem.topology.require_noneditable(func)[source]

Decorator: raise if topology is editable.

class mchem.topology.BondedAtoms(atoms: List[Atom])[source]

Bases: object

Ordered list of bonded atoms, with hash/equality invariant under reversal (for undirected paths).

__init__(atoms: List[Atom])[source]
Parameters:

atoms (List[Atom]) – Ordered sequence of bonded atoms.

extend(atoms)[source]
class mchem.topology.Topology(name: str, maxConnect: int = 5)[source]

Bases: object

Molecular topology: list of residues, bonds, and connectivity up to _maxConnect bonds.

__init__(name: str, maxConnect: int = 5)[source]
Parameters:
  • name (str) – Topology name (e.g. from PDB).

  • maxConnect (int, optional) – Maximum bond connectivity distance for generating bonded-atom lists.

property editable
setEditable()[source]
assignIndex()[source]
pruneBonds()[source]
property residues: List[Residue]
property numResidues
getResidueWithIdx(idx: int) Residue[source]
addResidue(*args, **kwargs)[source]
insertResidue(*args, **kwargs)[source]
removeResidue(*args, **kwargs)[source]
property numAtoms
property numBonds
atoms() Atom[source]
getAtomWithIdx(idx: int) Atom[source]
addAtomToResidue(*args, **kwargs)[source]
property bonds
addBond(*args, **kwargs)[source]
property connTable
property bondedAtoms
generateTopologicalInfo()[source]
matchTemplates()[source]
processDisulfideBond()[source]
property coordinates: ndarray

Solvation: add water and ions to a solute in a periodic box.

Uses Packmol as the external placement engine and computes box dimensions following OpenMM’s openmm.app.Modeller.addSolvent() conventions.

mchem.solvate.BoxShape

Allowed periodic box shapes.

alias of Literal[‘cube’, ‘dodecahedron’, ‘octahedron’]

mchem.solvate.solvate(topology: Topology, positions: ndarray, *, box_shape: Literal['cube', 'dodecahedron', 'octahedron'] = 'cube', buffer: float = 10.0, neutralize: bool = True, ionic_strength: float = 0.0, positive_ion: str = 'Na+', negative_ion: str = 'Cl-') Tuple[Topology, ndarray, Tuple[ndarray, ndarray, ndarray]][source]

Solvate a solute with water and optional ions in a periodic box.

Parameters:
  • topology (Topology) – Solute topology (no solvent).

  • positions (np.ndarray) – Atom positions, shape (N, 3), in Angstroms.

  • box_shape (BoxShape, optional) – Periodic box shape. Default is "cube".

  • buffer (float, optional) – Minimum padding in Angstroms between solute and box edge. Default is 10.0.

  • neutralize (bool, optional) – If True, add counterions to neutralize the system. Default is True.

  • ionic_strength (float, optional) – Target ionic strength in mol/L. Default is 0.0.

  • positive_ion (str, optional) – Positive ion type (e.g. "Na+"). Default is "Na+".

  • negative_ion (str, optional) – Negative ion type (e.g. "Cl-"). Default is "Cl-".

Returns:

(topology, positions, box_vectors) where box_vectors is a tuple of three np.ndarray vectors.

Return type:

tuple

Raises:

Physical constants and unit conversion factors for energy and length.

Periodic table data: element symbol, mass, and name.

class mchem.element.Element(atomicNum: int, symbol: str, mass: float, name: str)[source]

Bases: object

Dataclass for a chemical element.

atomicNum

Atomic number.

Type:

int

symbol

Element symbol (e.g. 'C', 'H').

Type:

str

mass

Atomic mass in amu.

Type:

float

name

Full element name.

Type:

str

atomicNum: int
symbol: str
mass: float
name: str
__init__(atomicNum: int, symbol: str, mass: float, name: str) None

Force field

Force field XML loading and term generation (see ForceField, Generator).

Force field definition from XML: atom types, residues, and force generators.

mchem.forcefield.base.str2float(string)[source]

Convert string to float.

mchem.forcefield.base.str2int(string)[source]

Convert string to int.

mchem.forcefield.base.float2str(number)[source]

Format number for XML (float or scientific).

mchem.forcefield.base.str2bool(string)[source]

Convert string to bool (True/False, case-insensitive).

mchem.forcefield.base.xmlele2str(xmlele: Element)[source]

Serialize XML element to string (unicode, stripped).

class mchem.forcefield.base.AtomType(name: str, atomClass: str)[source]

Bases: object

Force-field atom type with name and class.

name: str
atomClass: str
__init__(name: str, atomClass: str) None
class mchem.forcefield.base.ForceField(*files)[source]

Bases: object

Load and apply a force field from one or more XML files.

Parses AtomTypes, Residues, and force elements (e.g. AmoebaBondForce); uses Parsers to create generators that produce terms for a Topology.

__init__(*files)[source]

Load force field from one or more XML files.

Parameters:

*files (str or os.PathLike) – XML paths; resolved via processFileNames() if not found.

atomTypes: Dict[str, AtomType]
atomClasses: Dict[str, List[AtomType]]
processFileNames(files)[source]

Resolve file paths; search next to this module if not found. Returns list of paths.

property generators

List of force generators (one per force type from XML).

addGenerator(generator)[source]

Append a force generator.

getGeneratorWithClass(generatorClass)[source]

Return the first generator that is an instance of generatorClass, or None.

addGeneratorWithClass(generatorClass)[source]

Get or create a generator of the given class and return it.

addAtomType(atomTypeElement: Element)[source]

Register an atom type from an XML Type element (name, class).

loadAtomTypes()[source]

Load all AtomTypes from parsed XML trees.

loadAtomTypeDefs()[source]

Map residue atom names to atom types from Residues section (templates must exist).

loadForces()[source]

Parse force elements using Parsers and register generators.

assignAtomTypes(topology: Topology)[source]

Set atom type and class on each atom from residue templates.

findAtomTypes(element: Element, numAtoms: int)[source]

Resolve type/class attributes to list of atom type name tuples (supports wildcards).

createSystem(topology: Topology, **kwargs)[source]

Build a System with particles and all generator terms for the given topology.

getParameters(asJaxNumpy: bool = True)[source]
updateParameters(param: Dict[str, Any])[source]

Update each generator’s parameters from a dict keyed by generator class name.

exportAtomTypes()[source]

Return XML string for AtomTypes section.

exportAtomTypeDefs()[source]

Return XML string for Residues section (atom type defs).

save(path: PathLike)[source]

Write force field XML to file (atom types, residue defs, all force sections).

class mchem.forcefield.base.Generator(ff: ForceField, paramFields: List[str] = [], raiseError: bool = True)[source]

Bases: object

Base class for a force generator: holds parameters keyed by atom type/SMIRKS and creates terms for a topology.

__init__(ff: ForceField, paramFields: List[str] = [], raiseError: bool = True)[source]
getMetadata(key: str)[source]
setMetadata(key: str, value: Any)[source]
property paramFields: List[str]
property numParameters: int
addParameterField(name: str)[source]
checkParameter(paramDict: Dict[str, Any])[source]
addParameterWithAtomTypes(typeOrtypes, paramDict: Dict[str, Any])[source]
addParameterWithSmirks(smirks: str, paramDict: Dict[str, Any])[source]
setParameterWithIdx(paramDict: Dict[str, Any], paramIdx: int)[source]
getParameterIdxWithAtomType(typeQuery)[source]
getParameterWithAtomType(typeQuery)[source]
getParameterIdxWithSmirks(smirksQuery: str)[source]
getParameterWithSmirks(smirksQuery: str)[source]
getParameterWithIdx(idx: int)[source]
getParameters(asJaxNumpy: bool = False)[source]
updateParameters(param: Dict[str, Any])[source]
createTerms(topology: Topology, **kwargs)[source]

Build term list(s) for the given topology. Override in subclasses. Returns TermList or tuple of TermLists.

raise_exception(msg: str, raiseError: bool = True)[source]
exportParameterToStr()[source]

Force generators: create bonded/nonbonded terms from XML and topology.

class mchem.forcefield.generators.AmoebaBondGenerator(ff)[source]

Bases: Generator

Generator for AMOEBA bond terms (quartic) from AmoebaBondForce XML.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
addBond(bondElement: Element)[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.AmoebaAngleGenerator(ff)[source]

Bases: Generator

Generator for AMOEBA angle and in-plane angle terms from AmoebaAngleForce XML.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
addAngle(angleElement: Element)[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.AmoebaUreyBradleyGenerator(ff)[source]

Bases: Generator

Generator for AMOEBA Urey-Bradley terms from AmoebaUreyBradleyForce XML.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
addUreyBrad(ubElement: Element)[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.MultipoleGenerator(ff)[source]

Bases: Generator

Generator for atomic multipoles from AmoebaMultipoleForce / MBUCBMultipoleForce XML.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
static setAxisType(kz: int, kx: int, ky: int)[source]
addMultipole(mpoleElement: Element)[source]
exportParameterToStr()[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.IsotropicPolarizationGenerator(ff)[source]

Bases: Generator

Generator for isotropic polarizability terms from Polarize elements.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
addPolarize(polarElement: Element)[source]
setPolarizationGroup(topology: Topology)[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.AmoebaVdwGenerator(ff)[source]

Bases: Generator

Generator for AMOEBA buffered 14-7 VdW terms from AmoebaVdwForce XML.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
exportParameterToStr()[source]
addVdw(vdwElement: Element)[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.AmoebaStretchBendGenerator(ff)[source]

Bases: Generator

Generator for AMOEBA stretch-bend coupling from AmoebaStretchBendForce XML.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.AmoebaOutOfPlaneBendGenerator(ff)[source]

Bases: Generator

Generator for AMOEBA out-of-plane bend terms from AmoebaOutOfPlaneBendForce XML.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
addOutofPlaneBend(ele: Element)[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.AmoebaPiTorsionGenerator(ff)[source]

Bases: Generator

Generator for AMOEBA pi-torsion terms from AmoebaPiTorsionForce XML.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
addPiTorsion(ele: Element)[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.PeriodicTorsionGenerator(ff)[source]

Bases: Generator

Generator for periodic proper torsions from PeriodicTorsionForce XML.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
addProperTorsion(ele: Element)[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.AmoebaTorsionTorsionGenerator(ff)[source]

Bases: Generator

Generator for AMOEBA torsion-torsion (5-atom) terms and grids from AmoebaTorsionTorsionForce XML.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.AnisotropicPolarizationGenerator(ff)[source]

Bases: Generator

Generator for anisotropic polarizability from Polarize elements (MBUCB).

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
addPolarize(polarElement: Element)[source]
exportParameterToStr()[source]
setPolarizationGroup(topology: Topology)[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.MBUCBChargePenetrationGenerator(ff)[source]

Bases: Generator

Generator for MBUCB charge-penetration terms from ChargePenetration elements.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
addTerm(element: Element)[source]
exportParameterToStr()[source]
createTerms(topology: Topology, **kwargs)[source]
class mchem.forcefield.generators.MBUCBChargeTransferGenerator(ff)[source]

Bases: Generator

Generator for MBUCB charge-transfer terms from ChargeTransfer elements.

__init__(ff)[source]
static parseElement(element: Element, ff: ForceField)[source]
addTerm(element: Element)[source]
exportParameterToStr()[source]
createTerms(topology: Topology, **kwargs)[source]

Terms

Force term dataclasses and TermList (bonded, nonbonded).

Base container for force-term lists keyed by term class.

class mchem.terms.base.TermList(cls)[source]

Bases: list

List of force terms restricted to a single dataclass type.

Used by mchem.system.System to store terms per table (e.g. AmoebaBond).

__init__(cls)[source]
Parameters:

cls (type) – Dataclass type; only instances of this type can be appended.

property cls

Term dataclass type for this list.

append(item)[source]

Append a term; must be an instance of cls.

Bonded force terms: bonds, angles, Urey-Bradley, torsions, CMAP, AMOEBA-specific.

class mchem.terms.bonded.HarmonicBond(p0: int, p1: int, b0: float, kb: float, paramIdx: int = -1)[source]

Bases: object

Class for a bond with normal harmonic potential

p0: int
p1: int
b0: float
kb: float
paramIdx: int = -1
__init__(p0: int, p1: int, b0: float, kb: float, paramIdx: int = -1) None
class mchem.terms.bonded.AmoebaBond(p0: int, p1: int, b0: float, kb: float, cubic: float, quartic: float, paramIdx: int = -1)[source]

Bases: object

Class for a bond in AMOEBA force field, i.e. up to quartic polynomials

p0: int
p1: int
b0: float
kb: float
cubic: float
quartic: float
paramIdx: int = -1
__init__(p0: int, p1: int, b0: float, kb: float, cubic: float, quartic: float, paramIdx: int = -1) None
class mchem.terms.bonded.HarmonicAngle(p0: int, p1: int, p2: int, th0: float, kth: float, paramIdx: int = -1)[source]

Bases: object

Class for an angle with harmonic potential

p0: int
p1: int
p2: int
th0: float
kth: float
paramIdx: int = -1
__init__(p0: int, p1: int, p2: int, th0: float, kth: float, paramIdx: int = -1) None
class mchem.terms.bonded.AmoebaAngle(p0: int, p1: int, p2: int, th0: float, kth: float, cubic: float, quartic: float, pentic: float, sextic: float, paramIdx: int = -1)[source]

Bases: object

Class for an angle in AMOEBA force field (up to sextic polynomials)

p0: int
p1: int
p2: int
th0: float
kth: float
cubic: float
quartic: float
pentic: float
sextic: float
paramIdx: int = -1
__init__(p0: int, p1: int, p2: int, th0: float, kth: float, cubic: float, quartic: float, pentic: float, sextic: float, paramIdx: int = -1) None
class mchem.terms.bonded.AmoebaAngleInPlane(p0: int, p1: int, p2: int, p3: int, th0: float, kth: float, cubic: float, quartic: float, pentic: float, sextic: float, paramIdx: int = -1)[source]

Bases: object

Class for an in-plane angle in AMOEBA force field (up to sextic polynomials)

p0: int
p1: int
p2: int
p3: int
th0: float
kth: float
cubic: float
quartic: float
pentic: float
sextic: float
paramIdx: int = -1
__init__(p0: int, p1: int, p2: int, p3: int, th0: float, kth: float, cubic: float, quartic: float, pentic: float, sextic: float, paramIdx: int = -1) None
class mchem.terms.bonded.AmoebaStretchBend(p0: int, p1: int, p2: int, th0: float, b01: float, b02: float, kb1: float, kb2: float, paramIdx: int = -1)[source]

Bases: object

Class for stretch-bend coupling term in AMOEBA force field

p0: int
p1: int
p2: int
th0: float
b01: float
b02: float
kb1: float
kb2: float
paramIdx: int = -1
__init__(p0: int, p1: int, p2: int, th0: float, b01: float, b02: float, kb1: float, kb2: float, paramIdx: int = -1) None
class mchem.terms.bonded.AmoebaUreyBradley(p0: int, p1: int, p2: int, r0: float, fc: float, paramIdx: int = -1)[source]

Bases: object

Class for Urey-Bradley term used in AMOEBA force field

p0: int
p1: int
p2: int
r0: float
fc: float
paramIdx: int = -1
__init__(p0: int, p1: int, p2: int, r0: float, fc: float, paramIdx: int = -1) None
class mchem.terms.bonded.AmoebaOutOfPlaneBend(p0: int, p1: int, p2: int, p3: int, fc: float, paramIdx: int = -1)[source]

Bases: object

Class for out-of-plane bending term in AMOEBA force field

p0: int
p1: int
p2: int
p3: int
fc: float
paramIdx: int = -1
__init__(p0: int, p1: int, p2: int, p3: int, fc: float, paramIdx: int = -1) None
class mchem.terms.bonded.PeriodicTorsion(p0: int, p1: int, p2: int, p3: int, phase1: float, phase2: float, phase3: float, phase4: float, phase5: float, phase6: float, k1: float, k2: float, k3: float, k4: float, k5: float, k6: float, paramIdx: int = -1)[source]

Bases: object

Periodic torsion (dihedral) term with up to six Fourier components.

p0: int
p1: int
p2: int
p3: int
phase1: float
phase2: float
phase3: float
phase4: float
phase5: float
phase6: float
k1: float
k2: float
k3: float
k4: float
k5: float
k6: float
paramIdx: int = -1
__init__(p0: int, p1: int, p2: int, p3: int, phase1: float, phase2: float, phase3: float, phase4: float, phase5: float, phase6: float, k1: float, k2: float, k3: float, k4: float, k5: float, k6: float, paramIdx: int = -1) None
class mchem.terms.bonded.AmoebaStretchTorsion(p0: int, p1: int, p2: int, p3: int, k11: float, k12: float, k13: float, k21: float, k22: float, k23: float, k31: float, k32: float, k33: float, b01: float, b02: float, b03: float, phi01: float, phi02: float, phi03: float)[source]

Bases: object

Class for AMOEBA stretch-torsion coupling term

p0: int
p1: int
p2: int
p3: int
k11: float
k12: float
k13: float
k21: float
k22: float
k23: float
k31: float
k32: float
k33: float
b01: float
b02: float
b03: float
phi01: float
phi02: float
phi03: float
__init__(p0: int, p1: int, p2: int, p3: int, k11: float, k12: float, k13: float, k21: float, k22: float, k23: float, k31: float, k32: float, k33: float, b01: float, b02: float, b03: float, phi01: float, phi02: float, phi03: float) None
class mchem.terms.bonded.AmoebaAngleTorsion(p0: int, p1: int, p2: int, p3: int, k11: float, k12: float, k13: float, k21: float, k22: float, k23: float, th01: float, th02: float, phi01: float, phi02: float, phi03: float, paramIdx: int = -1)[source]

Bases: object

Class for AMOEBA angle-torsion coupling term

p0: int
p1: int
p2: int
p3: int
k11: float
k12: float
k13: float
k21: float
k22: float
k23: float
th01: float
th02: float
phi01: float
phi02: float
phi03: float
paramIdx: int = -1
__init__(p0: int, p1: int, p2: int, p3: int, k11: float, k12: float, k13: float, k21: float, k22: float, k23: float, th01: float, th02: float, phi01: float, phi02: float, phi03: float, paramIdx: int = -1) None
class mchem.terms.bonded.AmoebaPiTorsion(p0: int, p1: int, p2: int, p3: int, p4: int, p5: int, k: float, paramIdx: int = -1)[source]

Bases: object

Class for pi-torsion term in Amoeba

p0: int
p1: int
p2: int
p3: int
p4: int
p5: int
k: float
paramIdx: int = -1
__init__(p0: int, p1: int, p2: int, p3: int, p4: int, p5: int, k: float, paramIdx: int = -1) None
class mchem.terms.bonded.CMAPTable(phi: float, psi: float, ene: float, paramIdx: int = -1)[source]

Bases: object

Class for a CMAP tabulated torsion-torsion coupling term

phi: float
psi: float
ene: float
paramIdx: int = -1
__init__(phi: float, psi: float, ene: float, paramIdx: int = -1) None
class mchem.terms.bonded.CMAP(cmap: int, p0: int, p1: int, p2: int, p3: int, p4: int, p5: int, p6: int, p7: int)[source]

Bases: object

Class for CMAP torsion-torsion coupling term

cmap: int
p0: int
p1: int
p2: int
p3: int
p4: int
p5: int
p6: int
p7: int
__init__(cmap: int, p0: int, p1: int, p2: int, p3: int, p4: int, p5: int, p6: int, p7: int) None
class mchem.terms.bonded.AmoebaTorsionTorsionGrid(angle1: float, angle2: float, f: float, gridIdx: int)[source]

Bases: object

Class for an AMOEBA torsion-torsion 2D energy grid point

angle1: float
angle2: float
f: float
gridIdx: int
__init__(angle1: float, angle2: float, f: float, gridIdx: int) None
class mchem.terms.bonded.AmoebaTorsionTorsion(p0: int, p1: int, p2: int, p3: int, p4: int, gridIdx: int, nx: int, ny: int)[source]

Bases: object

Class for an AMOEBA torsion-torsion coupling term (5-atom chain)

p0: int
p1: int
p2: int
p3: int
p4: int
gridIdx: int
nx: int
ny: int
__init__(p0: int, p1: int, p2: int, p3: int, p4: int, gridIdx: int, nx: int, ny: int) None

Non-bonded terms: particles, VdW, multipoles, polarization, MBUCB.

class mchem.terms.nonbonded.MultipoleAxisType(value)[source]

Bases: Enum

Axis convention for multipole frame (ZThenX, Bisector, ZOnly, etc.).

ZThenX = 0
Bisector = 1
ZBisect = 2
ThreeFold = 3
ZOnly = 4
NoAxisType = 5
LastAxisTypeIndex = 6
class mchem.terms.nonbonded.Particle(idx: int, name: str, element: str, mass: float, resnum: int, resname: str, xx: float, xy: float, xz: float, vx: float = 0.0, vy: float = 0.0, vz: float = 0.0)[source]

Bases: object

Single particle: index, name, element, mass, residue info, position, optional velocity.

idx: int
name: str
element: str
mass: float
resnum: int
resname: str
xx: float
xy: float
xz: float
vx: float = 0.0
vy: float = 0.0
vz: float = 0.0
__init__(idx: int, name: str, element: str, mass: float, resnum: int, resname: str, xx: float, xy: float, xz: float, vx: float = 0.0, vy: float = 0.0, vz: float = 0.0) None
class mchem.terms.nonbonded.AmoebaVdw147(idx: int, epsilon: float, sigma: float, parentIdx: int = -1, reduction: float = 1.0, paramIdx: int = -1)[source]

Bases: object

AMOEBA buffered 14-7 VdW parameters (epsilon, sigma, reduction, optional parent).

idx: int
epsilon: float
sigma: float
parentIdx: int = -1
reduction: float = 1.0
paramIdx: int = -1
__init__(idx: int, epsilon: float, sigma: float, parentIdx: int = -1, reduction: float = 1.0, paramIdx: int = -1) None
class mchem.terms.nonbonded.Multipole(idx: int, c0: float, dx: float, dy: float, dz: float, qxx: float, qxy: float, qxz: float, qyy: float, qyz: float, qzz: float, axistype: int, kz: int = -1, kx: int = -1, ky: int = -1, paramIdx: int = -1)[source]

Bases: object

Atomic multipole (charge, dipole, quadrupole) and axis type / frame indices.

idx: int
c0: float
dx: float
dy: float
dz: float
qxx: float
qxy: float
qxz: float
qyy: float
qyz: float
qzz: float
axistype: int
kz: int = -1
kx: int = -1
ky: int = -1
paramIdx: int = -1
__init__(idx: int, c0: float, dx: float, dy: float, dz: float, qxx: float, qxy: float, qxz: float, qyy: float, qyz: float, qzz: float, axistype: int, kz: int = -1, kx: int = -1, ky: int = -1, paramIdx: int = -1) None
class mchem.terms.nonbonded.IsotropicPolarization(idx: int, alpha: float, thole: float, grp: list, paramIdx: int = -1)[source]

Bases: object

Isotropic polarizability (alpha, thole) and polarization group indices.

idx: int
alpha: float
thole: float
grp: list
paramIdx: int = -1
__init__(idx: int, alpha: float, thole: float, grp: list, paramIdx: int = -1) None
class mchem.terms.nonbonded.AnisotropicPolarization(idx: int, alphaxx: float, alphaxy: float, alphaxz: float, alphayy: float, alphayz: float, alphazz: float, thole: float, grp: list = <factory>, paramIdx: int = -1)[source]

Bases: object

Anisotropic polarizability tensor (alpha_xx, alpha_xy, …) and thole/group.

idx: int
alphaxx: float
alphaxy: float
alphaxz: float
alphayy: float
alphayz: float
alphazz: float
thole: float
grp: list
paramIdx: int = -1
__init__(idx: int, alphaxx: float, alphaxy: float, alphaxz: float, alphayy: float, alphayz: float, alphazz: float, thole: float, grp: list = <factory>, paramIdx: int = -1) None
class mchem.terms.nonbonded.MBUCBChargePenetration(idx: int, z: float, alpha: float, beta: float, paramIdx: int = -1)[source]

Bases: object

MBUCB charge-penetration correction parameters (z, alpha, beta).

idx: int
z: float
alpha: float
beta: float
paramIdx: int = -1
__init__(idx: int, z: float, alpha: float, beta: float, paramIdx: int = -1) None
class mchem.terms.nonbonded.MBUCBChargeTransfer(idx: int, d: float, b: float, alpha: float, paramIdx: int = -1)[source]

Bases: object

MBUCB charge-transfer parameters (d, b, alpha).

idx: int
d: float
b: float
alpha: float
paramIdx: int = -1
__init__(idx: int, d: float, b: float, alpha: float, paramIdx: int = -1) None
class mchem.terms.nonbonded.PairList(p0: int, p1: int)[source]

Bases: object

Class for atom pairs that has to specially treated in non-bonded calculations

p0: int
p1: int
energy() float[source]
__init__(p0: int, p1: int) None

File formats

File format loaders (e.g. PDB) that produce mchem.topology.Topology.

PDB file I/O: read/write Topology from/to ATOM/HETATM records.

mchem.fileformats.pdb.read_pdb_box(fname: PathLike) Tuple[float, float, float, float, float, float] | None[source]

Read periodic box dimensions from the first CRYST1 record in a PDB file.

Parameters:

fname (os.PathLike) – Path to the PDB file.

Returns:

Cell lengths a, b, c in Angstroms and angles alpha, beta, gamma in degrees. None if no CRYST1 record is present.

Return type:

tuple of (a, b, c, alpha, beta, gamma) or None

mchem.fileformats.pdb.load_pdb(fname: PathLike) Topology[source]

Load a topology from a PDB file (ATOM/HETATM), match residue templates, and return editable topology.

Parameters:

fname (os.PathLike) – Path to the PDB file.

Returns:

Topology with residues and atoms; templates are matched so bonds and standard names are set.

Return type:

Topology

mchem.fileformats.pdb.write_pdb(fname: PathLike, topology: Topology, positions: ndarray, box_vectors: Tuple[ndarray, ndarray, ndarray] | None = None) None[source]

Write a PDB file from a topology and positions.

Parameters:
  • fname (os.PathLike) – Output PDB file path.

  • topology (Topology) – Molecular topology.

  • positions (np.ndarray) – Atom positions, shape (N, 3), in Angstroms.

  • box_vectors (tuple of np.ndarray, optional) – Three box vectors (v1, v2, v3). If provided, a CRYST1 record is written at the top of the file.