pysimm.system module

class pysimm.system.Angle(**kwargs)[source]

Bases: pysimm.utils.Item

Angle between particles a, b, and c

a–b–c

Objects inheriting from Item can contain arbitrary data. Keyword arguments are assigned as attributes. Attributes usually used are given below.

a

Particle object involved in angle

b

Particle object involved in angle (middle particle)

c

Particle object involved in angle

type

AngleType object reference

angle(radians=False)[source]

Calculate angle.

Parameters:radians – True to return value in radians (default: False)
Returns:Angle between Particle a, b, and c
class pysimm.system.AngleType(**kwargs)[source]

Bases: pysimm.utils.Item

Objects inheriting from Item can contain arbitrary data. Keyword arguments are assigned as attributes. Attributes usually used are given below.

k

harmonic angle bend force constant (kcal/mol/radian^2)

theta0

angle equilibrium value (degrees)

name

force field angle type name

form(style='harmonic', d_range=None)[source]

Returns data to plot functional form for the potential energy with the given style.

Parameters:style – string for pair style of AngleType (harmonic, class2, charmm)
Returns:x, y for plotting functional form (energy vs angle)
classmethod guess_style(nparam)[source]
classmethod parse_lammps(line, style)[source]
write_lammps(style='harmonic', cross_term=None)[source]

Formats a string to define angle type coefficients for a LAMMPS data file given the provided style.

Parameters:
  • style – string for pair style of AngleType (harmonic, class2, charmm)
  • cross_term – type of class2 cross term to write (default=None) - BondBond - BondAngle
Returns:

LAMMPS formatted string with angle coefficients

class pysimm.system.Bond(**kwargs)[source]

Bases: pysimm.utils.Item

Bond between particle a and b

a–b

Objects inheriting from Item can contain arbitrary data. Keyword arguments are assigned as attributes. Attributes usually used are given below.

a

Particle object involved in bond

b

Particle object involved in bond

type

BondType object reference

distance()[source]

Calculates distance between Particle a and Particle b in this Bond object. Sets distance to dist attribute of self. Does not consider pbc.

Parameters:None
Returns:Distance between Particle a and Particle b (not considering pbc)
get_other_particle(p)[source]
class pysimm.system.BondType(**kwargs)[source]

Bases: pysimm.utils.Item

Objects inheriting from Item can contain arbitrary data. Keyword arguments are assigned as attributes. Attributes usually used are given below.

k

harmonic bond force constant (kcal/mol/A^2)

r0

bond equilibrium distance (Angstrom)

name

force field bond type name

form(style='harmonic', d_range=None)[source]

Returns data to plot functional form for the potential energy with the given style.

Parameters:style – string for pair style of BondType (harmonic, class2)
Returns:x, y for plotting functional form (energy vs distance)
classmethod guess_style(nparam)[source]
classmethod parse_lammps(line, style)[source]
write_lammps(style='harmonic')[source]

Formats a string to define bond type coefficients for a LAMMPS data file given the provided style.

Parameters:style – string for pair style of BondType (harmonic, class2)
Returns:LAMMPS formatted string with bond coefficients
class pysimm.system.Dihedral(**kwargs)[source]

Bases: pysimm.utils.Item

Dihedral between particles a, b, c, and d

a–b–c–d

Objects inheriting from Item can contain arbitrary data. Keyword arguments are assigned as attributes. Attributes usually used are given below.

a

Particle object involved in dihedral

b

Particle object involved in dihedral (middle particle)

c

Particle object involved in dihedral (middle particle)

d

Particle object involved in dihedral

type

DihedralType object reference

class pysimm.system.DihedralType(**kwargs)[source]

Bases: pysimm.utils.Item

Objects inheriting from Item can contain arbitrary data. Keyword arguments are assigned as attributes. Attributes usually used are given below.

k

dihedral energy barrier (kcal/mol)

d

minimum (+1 or -1)

n

multiplicity (integer >= 0)

name

force field dihedral type name

form(style='harmonic', d_range=None)[source]

Returns data to plot functional form for the potential energy with the given style.

Parameters:style – string for pair style of DihedralType (harmonic, class2, fourier)
Returns:x, y for plotting functional form (energy vs angle)
classmethod guess_style(nparam)[source]
classmethod parse_lammps(line, style)[source]
write_lammps(style='harmonic', cross_term=None)[source]

Formats a string to define dihedral type coefficients for a LAMMPS data file given the provided style.

Parameters:
  • style – string for pair style of DihedralType (harmonic, class2, fourier)
  • cross_term – type of class2 cross term to write (default=None) - MiddleBond - EndBond - Angle - AngleAngle - BondBond13
Returns:

LAMMPS formatted string with dihedral coefficients

class pysimm.system.Dimension(**kwargs)[source]

Bases: pysimm.utils.Item

Objects inheriting from Item can contain arbitrary data. Keyword arguments are assigned as attributes. Attributes usually used are given below.

xlo

minimum value in x dimension

xhi

maximum value in x dimension

ylo

minimum value in y dimension

yhi

maximum value in y dimension

zlo

minimum value in z dimension

zhi

maximum value in z dimension

dx

distance in x dimension

dy

distance in y dimension

dz

distance in z dimension

check()[source]
dx
dy
dz
size()[source]
translate(x, y, z)[source]

Shifts box bounds by x, y, z.

Parameters:
  • x – distance to shift box bounds in x direction
  • y – distance to shift box bounds in y direction
  • z – distance to shift box bounds in z direction
Returns:

None

class pysimm.system.Improper(**kwargs)[source]

Bases: pysimm.utils.Item

Improper dihedral around particle a, bonded to b, c, and d

b
a–d
c

Objects inheriting from Item can contain arbitrary data. Keyword arguments are assigned as attributes. Attributes usually used are given below.

a

Particle object involved in improper (middle particle)

b

Particle object involved in improper

c

Particle object involved in improper

d

Particle object involved in improper

type

ImproperType object reference

class pysimm.system.ImproperType(**kwargs)[source]

Bases: pysimm.utils.Item

Objects inheriting from Item can contain arbitrary data. Keyword arguments are assigned as attributes. Attributes usually used are given below.

k

improper energy barrier (kcal/mol)

x0

equilibrium value (degrees)

name

force field improper type name

form(style='harmonic', d_range=None)[source]

Returns data to plot functional form for the potential energy with the given style.

Parameters:style – string for pair style of ImproperType (harmonic, cvff)
Returns:x, y for plotting functional form (energy vs angle)
classmethod guess_style(nparam)[source]
classmethod parse_lammps(line, style)[source]
write_lammps(style='harmonic', cross_term=None)[source]

Formats a string to define improper type coefficients for a LAMMPS data file given the provided style.

Parameters:
  • style – string for pair style of ImproperType (harmonic, class2, cvff)
  • cross_term – type of class2 cross term to write (default=None) - AngleAngle
Returns:

LAMMPS formatted string with dihedral coefficients

class pysimm.system.Molecule(**kwargs)[source]

Bases: pysimm.system.System

Very similar to System, but requires less information

class pysimm.system.Particle(**kwargs)[source]

Bases: pysimm.utils.Item

Objects inheriting from Item can contain arbitrary data. Keyword arguments are assigned as attributes. Attributes usually used are given below.

x

x coordinate

y

y coordinate

z

z coordinate

charge

partial charge

type

ParticleType object reference

check(style='full')[source]
coords()[source]
delete_bonding(s)[source]

Iterates through s.bonds, s.angles, s.dihedrals, and s.impropers and removes those which contain this Particle.

Parameters:sSystem object from which bonding objects will be removed
Returns:None
translate(dx, dy, dz)[source]

Shifts Particle position by dx, dy, dz.

Parameters:
  • dx – distance to shift in x direction
  • dy – distance to shift in y direction
  • dz – distance to shift in z direction
Returns:

None

class pysimm.system.ParticleType(**kwargs)[source]

Bases: pysimm.utils.Item

Objects inheriting from Item can contain arbitrary data. Keyword arguments are assigned as attributes. Attributes usually used are given below.

sigma

LJ sigma value (Angstrom)

epsilon

LJ epsilon value (kcal/mol)

elem

element abbreviation, i.e. ‘H’ for Hydrogen, ‘Cl’ for Chlorine

name

force field particle type name

form(style='lj_12-6', d_range=None)[source]

Returns data to plot functional form for the potential energy with the given style.

Parameters:style – string for pair style of ParticleType (lj_12-6, lj_9-6, buck)
Returns:x, y for plotting functional form (energy vs distance)
classmethod guess_style(nparam)[source]
classmethod parse_lammps(line, style)[source]
write_lammps(style='lj')[source]

Formats a string to define particle type coefficients for a LAMMPS data file given the provided style.

Parameters:style – string for pair style of ParticleType (lj, class2, mass, buck)
Returns:LAMMPS formatted string with pair coefficients
class pysimm.system.System(**kwargs)[source]

Bases: object

Object representation of molecular system. Contains information required for molecular simulation.

dim

Dimension object reference

particles

ItemContainer for Particle organization

particle_types

ItemContainer for ParticleType organization

bonds

ItemContainer for Bond organization

bond_types

ItemContainer for BondType organization

angles

ItemContainer for Angle organization

angle_types

ItemContainer for AngleType organization

dihedrals

ItemContainer for Dihedral organization

dihedral_types

ItemContainer for DihedralType organization

impropers

ItemContainer for Improper organization

improper_types

ItemContainer for ImproperType organization

molecules

ItemContainer for Molecule organization

add(other, **kwargs)[source]

Add other System to this. Optionally remove duplicate types (default behavior).

Parameters:
  • otherSystem object to add
  • unique_types (optional) – Remove duplicate types and reassign references to existing types (True)
  • change_dim (optional) – Update Dimension object so that Particle objects do not exist outside of Dimension extremes (True)
  • update_properties (optional) – Update system-wide mass, volume, density, center of gravity, and velocity properties (True)
add_angle(a=None, b=None, c=None, f=None)[source]

Add Angle to system between three particles

Parameters:
  • aParticle involved in new Angle
  • bParticle involved in new Angle (middle particle)
  • cParticle involved in new Angle
  • fForcefield object from which new force field type will be retrieved
Returns:

None

add_bond(a=None, b=None, f=None)[source]

Add Bond to system between two particles

Parameters:
  • aParticle involved in new Bond
  • bParticle involved in new Bond
  • fForcefield object from which new force field type will be retrieved
Returns:

None

add_dihedral(a=None, b=None, c=None, d=None, f=None)[source]

Add Dihedral to system between four particles

Parameters:
Returns:

None

add_improper(a=None, b=None, c=None, d=None, f=None)[source]

Add Improper to system between four particles

Parameters:
  • aParticle involved in new Improper (middle particle)
  • bParticle involved in new Improper
  • cParticle involved in new Improper
  • dParticle involved in new Improper
  • fForcefield object from which new force field type will be retrieved
Returns:

None

add_particle(p)[source]

Add new Particle to System.

Parameters:p – new Particle object to be added
Returns:None
add_particle_bonded_to(p, p0, f=None, sep=1.5)[source]

Add new Particle to System bonded to p0 and automatically update new forcefield types

Parameters:
  • p – new Particle object to be added
  • p0 – original Particle object in System to which p will be bonded
  • fForcefield object from which new force field types will be retrieved
Returns:

new Particle being added to system for convenient reference

add_particle_bonding()[source]

Update Particle objects such that Particle.bonded_to contains other Particle objects invloved in bonding

Parameters:None
Returns:None
apply_charges(f, charges='default')[source]

Applies charges derived using method provided by user. Defaults to ‘default’. Calls assign_charges() method of forcefield object provided.

Parameters:
  • fForcefield object
  • charges – type of charges to be applied default=’default’
Returns:

None

apply_forcefield(f, charges='default', set_box=True, box_padding=10, update_ptypes=False, skip_ptypes=False)[source]

Applies force field data to System based on typing rules defined in Forcefield object f

Parameters:
  • fForcefield object from which new force field type will be retrieved
  • charges – type of charges to be applied default=’default’
  • set_box – Update simulation box information based on particle positions default=True
  • box_padding – Add padding to simulation box if updating dimensions default=10 (Angstroms)
  • update_ptypes – If True, update particle types based on current ParticleType names default=False
  • skip_ptypes – if True, do not change particle types
Returns:

None

center(what='particles', at=[0, 0, 0], move_both=True)[source]

Centers particles center of geometry or simulation box at given coordinate. A vector is defined based on the current coordinate for the center of either the particles or the simulation box and the “at” parameter. This shift vector is applied to the entity defined by the “what” parameter. Optionally, both the particles and the box can be shifted by the same vector.

Parameters:
  • what – what is being centered: “particles” or “box”
  • at – new coordinate for center of particles or box
  • move_both – if True, determines vector for shift defined by “what” and “at” parameters, and applies shift to both particles and box. If false, only shift what is defined by “what” parameter.
Returns:

None

center_at_origin()[source]

DEPRECATED: Use System.center(‘particles’, [0, 0, 0], True) instead

Parameters:None
Returns:None
center_system()[source]

DEPRECATED: Use System.center(‘box’, [0, 0, 0], True) instead

Parameters:None
Returns:None
check_forcefield()[source]

Iterates through particles and prints the following:

tag type name type element type description bonded elements

Parameters:None
Returns:None
check_items()[source]

Checks particles, bonds, angles, dihedrals, impropers, and molecules containers and raises exception if the length of items in the container does not equal the count property

Parameters:None
Returns:None
consolidate_types()[source]

Removes duplicate types and reassigns references

Parameters:None
Returns:None
copy(rotate_x=None, rotate_y=None, rotate_z=None, dx=0, dy=0, dz=0)[source]

Create duplicate System object. Default behavior does not modify particle positions.

Parameters:
  • rotate_x – rotate duplicate system around x axis by this value (radians)
  • rotate_y – rotate duplicate system around y axis by this value (radians)
  • rotate_z – rotate duplicate system around z axis by this value (radians)
  • dx – translate duplicate system in x dimension by this value (Angstrom)
  • dy – translate duplicate system in y dimension by this value (Angstrom)
  • dz – translate duplicate system in z dimension by this value (Angstrom)
distance(p1, p2)[source]

Calculate distance between two particles considering pbc.

Parameters:
Returns:

distance between particles considering pbc

make_linker_types()[source]

Identifies linker particles and creates duplicate Particle.linker attribute. New ParticleType name is prepended with [H or T]L@ to designate head or tail linker

Parameters:None
Returns:None
make_new_bonds(p1=None, p2=None, f=None, angles=True, dihedrals=True, impropers=True)[source]

Makes new bond between two particles and updates new force field types

Parameters:
  • p1Particle object involved in new bond
  • p2Particle object involved in new bond
  • fForcefield object from which new force field types will be retrieved
  • angles – True to update new angles default=True
  • dihedrals – True to update new dihedrals default=True
  • impropers – True to update new impropers default=True
Returns:

None

objectify()[source]

Set references for Bond, Angle, Dihedral, Improper objects. For example, if read from file that bond #1 is between particle 1 and 2 set Bond.a to particles[1], etc.

Parameters:None
Returns:None
particles_df(columns=['tag', 'x', 'y', 'z', 'q'], index='tag', extras=[])[source]
quality(tolerance=0.1)[source]

Attemps to assess quality of System based on bond lengths in unwrapped system.

Parameters:tolerance – fractional value of equilibrium bond length that is acceptable
Returns:number of bonds in system outside tolerance
read_lammps_dump(fname)[source]

Updates particle positions and box size from LAMMPS dump file. Assumes following format for each atom line:

tag charge xcoord ycoord zcoord xvelocity yvelocity zvelocity

Parameters:fname – LAMMPS dump file
Returns:None
read_lammpstrj(trj, frame=1)[source]

Updates particle positions and box size from LAMMPS trajectory file at given frame.

Assumes one of following formats for each atom line:

tag xcoord ycoord zcoord

OR

tag type_id xcoord ycoord zcoord

OR

tag type_id xcoord ycoord zcoord ximage yimage zimage

Parameters:
  • trj – LAMMPS trajectory file
  • frame – sequential frame number (not LAMMPS timestep) default=1
Returns:

None

read_type_names(types_file)[source]

Update ParticleType names from file.

Parameters:types_file – type dictionary file name
Returns:None
read_xyz(xyz, frame=1, quiet=False)[source]

Updates particle positions and box size from xyz file at given frame

Parameters:
  • xyz – xyz trajectory file
  • frame – sequential frame number default=1
  • quiet – True to print status default=False
Returns:

None

remove_linker_types()[source]

Reassigns Particle.type references to original ParticleType objects without linker prepend

Parameters:None
Returns:None
remove_spare_bonding(update_tags=True)[source]

Removes bonds, angles, dihedrals and impropers that reference particles not in System.particles

Parameters:update_tags – True to update all tags after removal of bonding items default=True
rotate(around=None, theta_x=0, theta_y=0, theta_z=0, rot_matrix=None)[source]

* REQUIRES NUMPY *

Rotates entire system around given Particle by user defined angles

Parameters:
  • aroundParticle around which System will be rotated default=None
  • theta_x – angle around which system will be rotated on x axis
  • theta_y – angle around which system will be rotated on y axis
  • theta_z – angle around which system will be rotated on z axis
  • rot_matrix – rotation matrix to use for rotation
Returns:

None

set_atomic_numbers()[source]

Updates ParticleType objects with atomic number based on ParticleType.elem

Parameters:None
Returns:None
set_box(padding=0.0, center=True)[source]

Update System.dim with user defined padding. Used to construct a simulation box if it doesn’t exist, or adjust the size of the simulation box following system modifications.

Parameters:
  • padding – add padding to all sides of box (Angstrom)
  • center – if True, place center of box at origin default=True
Returns:

None

set_charge()[source]

Sets total charge of all Particle objects in System.particles

Parameters:None
Returns:None
set_cog()[source]

Calculate center of gravity of System and assign to System.cog

Parameters:None
Returns:None
set_density()[source]

Calculate density of System from mass and volume

Parameters:None
Returns:None
set_excluded_particles(bonds=True, angles=True, dihedrals=True)[source]

Updates Particle object such that Particle.excluded_particles contains other Particle objects involved in 1-2, 1-3, and/or 1-4 interactions

Parameters:
  • bonds – exclude particles involved in 1-2 interactions
  • angles – exclude particles involved in 1-3 interactions
  • dihedrals – exclude particles involved in 1-4 interactions
set_frac_free_volume(v_void=None)[source]

Calculates fractional free volume from void volume and bulk density

Parameters:v_void – void volume if not defined in System.void_volume default=None
Returns:None
set_mass()[source]

Set total mass of particles in System

Parameters:None
Returns:None
set_mm_dist(molecules=None)[source]

Calculate molecular mass distribution (mainly for polymer systems). Sets System.mw, System.mn, and System.disperisty

Parameters:moleculesItemContainer of molecules to calculate distributions defaul=’all’
Returns:None
set_references()[source]

Set object references when System information read from text file. For example, if bond type value 2 is read from file, set Bond.type to bond_types[2]

Parameters:None
Returns:None
set_velocity()[source]

Calculate total velocity of particles in System

Parameters:None
Returns:None
set_volume()[source]

Set volume of System based on Dimension

Parameters:None
Returns:None
shift_particles(shiftx, shifty, shiftz)[source]

Shifts all particles by shiftx, shifty, shiftz. Recalculates cog.

Parameters:
  • shiftx – distance to shift particles in x direction
  • shifty – distance to shift particles in y direction
  • shiftz – distance to shift particles in z direction
Returns:

None

shift_to_origin()[source]

Shifts simulation box to begin at origin. i.e. xlo=ylo=zlo=0

Parameters:None
Returns:None
unite_atoms()[source]
unwrap()[source]

Unwraps Particle images such that no bonds cross box edges.

Parameters:None
Returns:None
update_ff_types_from_ac(ff, acname)[source]

Updates ParticleType objects in system using type names given in antechamber (ac) file. Retrieves type from System if possible, then searches force field provided by ff.

Parameters:
  • ff – forcefield to search for Type objects
  • acname – ac filename containing type names
Returns:

None

update_particle_types_from_forcefield(f)[source]

pysimm.system.System.update_types_from_forcefield

Updates ParticleType data from Forcefield object f based on ParticleType.name

Parameters:fForcefield object reference
Returns:None
update_tags()[source]

Update Item tags in ItemContainer objects to preserve continuous tags. Removes all objects and then reinserts them.

Args:
None
Returns:
None
update_types(ptypes, btypes, atypes, dtypes, itypes)[source]

Updates type objects from a given list of types.

Parameters:
  • ptypes – list of ParticleType objects from which to update
  • btypes – list of BondType objects from which to update
  • atypes – list of AngleType objects from which to update
  • dtypes – list of DihedralType objects from which to update
  • itypes – list of ImproperType objects from which to update
visualize(vis_exec='vmd', **kwargs)[source]

Visualize system in third party software with given executable. Software must accept pdb or xyz as first command line argument.

Parameters:
  • vis_exec – executable to launch visualization software default=’vmd’
  • unwrap (optional) – if True, unwrap System first default=None
  • format (optional) – set format default=’xyz’
Returns:

None

viz(**kwargs)[source]
wrap()[source]

Wrap Particle images into box defined by Dimension object. Ensure particles are contained within simulation box.

Parameters:None
Returns:None
write_chemdoodle_json(outfile, **kwargs)[source]

Write System data in chemdoodle json format

Parameters:outfile – where to write data, file name or ‘string’
Returns:None or string of data file if out_data=’string’
write_cssr(outfile='data.cssr', **kwargs)[source]

Write System data in cssr format file format: line, format, contents

1: 38X, 3F8.3 : - length of the three cell parameters (a, b, and c) in angstroms.

2: 21X, 3F8.3, 4X, ‘SPGR =’, I3, 1X, A11 : - a, b, g in degrees, space group number, space group name.

3: 2I4, 1X, A60 : - Number of atoms stored, coordinate system flag (0=fractional, 1=orthogonal coordinates in Angstrom), first title.

4: A53 : - A line of text that can be used to describe the file.

5: I4, 1X, A4, 2X, 3(F9.5.1X), 8I4, 1X, F7.3 : - Atom serial number, atom name, x, y, z coordinates, bonding connectivities (max 8), charge.

Note: The atom name is a concatenation of the element symbol and the atom serial number.

Keyword Arguments:
 
  • outfile – where to write data, file name or ‘string’ (default: ‘data.cssr’)
  • frac – flag for using fractional coordinates (default: False)
  • aname – flag for using element as an atom name (if it is True); else will use atom type name (default: False)
Returns:

None or string of data file if out_data=’string’

write_lammps(out_data, **kwargs)[source]

Write System data formatted for LAMMPS

Parameters:out_data – where to write data, file name or ‘string’
Returns:None or string if data file if out_data=’string’
write_lammps_mol(out_data)[source]

Write System data formatted as LAMMPS molecule template

Parameters:out_data – where to write data, file name or ‘string’
Returns:None or string if data file if out_data=’string’
write_mol(outfile='data.mol')[source]

Write System data in mol format

Parameters:outfile – where to write data, file name or ‘string’
Returns:None or string of data file if out_data=’string’
write_pdb(outfile='data.pdb', type_names=True)[source]

Write System data in pdb format

Parameters:outfile – where to write data, file name or ‘string’
Returns:None or string of data file if out_data=’string’
write_xyz(outfile='data.xyz', **kwargs)[source]

Write System data in xyz format

Parameters:outfile – where to write data, file name or ‘string’
Returns:None or string of data file if out_data=’string’
write_yaml(file_)[source]

Write System data in yaml format

Parameters:outfile – file name to write data
Returns:None
zero_charge()[source]

Enforces total System charge to be 0.0 by subtracting excess charge from last particle

Parameters:None
Returns:None
zero_velocity()[source]

Enforce zero shift velocity in System

Parameters:None
Returns:None
pysimm.system.compare(s1, s2)[source]
pysimm.system.distance_to_origin(p)[source]

Calculates distance of particle to origin.

Parameters:p – Particle object with x, y, and z attributes
Returns:Distance of particle to origin
pysimm.system.get_types(*arg, **kwargs)[source]

Get unique type names from list of systems

Parameters:write (optional) – if True, write types dictionary to filename
Returns:(ptypes, btypes, atypes, dtypes, itypes) * for use with update_types *
pysimm.system.random() → x in the interval [0, 1).
pysimm.system.read_ac(ac_file)[source]

Interprets ac file and creates System object

Parameters:ac_file – ac file name
Returns:System object
pysimm.system.read_chemdoodle_json(file_, **kwargs)[source]

pysimm.system.read_xyz

Interprets xyz file and creates System object

Parameters:
  • file – xyz file name
  • quiet (optional) – if False, print status
Returns:

System object

pysimm.system.read_cml(cml_file, **kwargs)[source]

Interprets cml file and creates System object

Parameters:
  • cml_file – cml file name
  • linkers (optional) – if True, use spinMultiplicity to determine linker default=None
Returns:

System object

pysimm.system.read_lammps(data_file, **kwargs)[source]

Interprets LAMMPS data file and creates System object

Parameters:
  • data_file – LAMMPS data file name
  • quiet (optional) – if False, print status
  • atom_style (optional) – option to let user override (understands charge, molecular, full)
  • pair_style (optional) – option to let user override
  • bond_style (optional) – option to let user override
  • angle_style (optional) – option to let user override
  • dihedral_style (optional) – option to let user override
  • improper_style (optional) – option to let user override
  • set_types (optional) – if True, objectify default=True
  • name (optional) – provide name for system
Returns:

System object

pysimm.system.read_mol(mol_file, type_with=None, version='V2000')[source]

Interprets mol file and creates System object

Parameters:
  • mol_file – mol file name
  • f (optional) – Forcefield object to get data from
  • version – version of mol file to expect default=’V2000’
Returns:

System object

pysimm.system.read_pdb(pdb_file)[source]

Interprets pdb file and creates System object

Parameters:pdb_file – pdb file name
Returns:System object
pysimm.system.read_prepc(prec_file)[source]

Interprets prepc file and creates System object

Parameters:prepc_file – ac file name
Returns:System object
pysimm.system.read_pubchem_cid(cid, type_with=None)[source]

pysimm.system.read_pubchem_smiles

Interface with pubchem restful API to create molecular system from SMILES format

Parameters:
  • smiles – smiles formatted string of molecule
  • type_withForcefield object to type with default=None
Returns:

System object

pysimm.system.read_pubchem_smiles(smiles, quiet=False, type_with=None)[source]

Interface with pubchem restful API to create molecular system from SMILES format

Parameters:
  • smiles – smiles formatted string of molecule
  • type_withForcefield object to type with default=None
Returns:

System object

pysimm.system.read_xyz(file_, **kwargs)[source]

Interprets xyz file and creates System object

Parameters:
  • file – xyz file name
  • quiet (optional) – if False, print status
Returns:

System object

pysimm.system.read_yaml(file_, **kwargs)[source]

Interprets yaml file and creates System object

Parameters:file – yaml file name
Returns:System object
pysimm.system.replicate(ref, nrep, s_=None, density=0.3, rand=True, print_insertions=True)[source]

Replicates list of System objects into new (or exisintg) System. Can be random insertion.

Parameters:
  • ref – reference :class:`~pysimm.system.System`(s) (this can be a list)
  • nrep – number of insertions to perform (can be list but must match length of ref)
  • sSystem into which insertions will be performed default=None
  • density – density of new System default=0.3 (set to None to not change box)
  • rand – if True, random insertion is performed
  • print_insertions – if True, update screen with number of insertions