API Reference¶
Classes¶
- class pymolar.AnalysisTask¶
Bases:
objectBase class for trajectory processing tasks.
Subclass and implement register_args, pre_process, process_frame, and post_process. The constructor parses CLI arguments, streams trajectory frames, and calls the hooks in processing order.
- __init__()¶
Run the full analysis pipeline from command-line arguments.
- register_args(parser)¶
Register task-specific CLI arguments on the provided parser.
- pre_process()¶
Hook called once before the first processed frame.
- process_frame()¶
Hook called for each processed frame.
- post_process()¶
Hook called once after all frames are processed.
- class pymolar.Atom¶
Bases:
objectMutable atom container.
Example
import pymolar a = pymolar.Atom() a.name = "CA" a.resname = "ALA" a.resid = 1 a.chain = 'A' a.mass = 12.011
- atomic_number¶
Atomic number.
- bfactor¶
Temperature factor (B-factor).
- chain¶
Chain identifier.
- charge¶
Atom charge.
- mass¶
Atomic mass.
- name¶
Atom name.
- occupancy¶
Occupancy value.
- resid¶
Residue identifier.
- resname¶
Residue name.
- type_id¶
Force-field atom type identifier.
- type_name¶
Force-field atom type name.
- class pymolar.FileHandler(fname, mode)¶
Bases:
objectReader/writer for topology and trajectory files.
Example
import pymolar as molar fh = molar.FileHandler("traj.xtc", "r") for st in fh: print(st.time)
- file_name¶
File path associated with this handler.
- Returns:
File path.
- Return type:
str
- read()¶
Read topology and the next state frame.
- skip_to_frame(fr)¶
Seek reader to frame index.
- Parameters:
fr – Target frame index.
- Returns:
None.- Return type:
None
- skip_to_last()¶
Seek the reader to the last frame.
- Returns:
None.- Return type:
None
- skip_to_time(t)¶
Seek reader to simulation time.
- Parameters:
t – Target simulation time.
- Returns:
None.- Return type:
None
- write(data)¶
Write a
System,Sel, or(Topology, State)tuple.- Parameters:
data – Data object to write.
- Returns:
None.- Return type:
None
- write_state(data)¶
Write state from
System,Sel, orState.- Parameters:
data – Source object containing state.
- Returns:
None.- Return type:
None
- write_topology(data)¶
Write topology from
System,Sel, orTopology.- Parameters:
data – Source object containing topology.
- Returns:
None.- Return type:
None
- class pymolar.FileStats¶
Bases:
objectRuntime IO statistics collected by
FileHandler.Example
import pymolar fh = pymolar.FileHandler("traj.xtc", "r") stats = fh.stats print(stats.frames_processed, stats.cur_t)
- cur_t¶
Time value of the current frame.
- Returns:
Current frame time.
- Return type:
float
- elapsed_time¶
Total wall-clock time spent in IO.
- Returns:
Elapsed duration.
- Return type:
datetime.timedelta
- frames_processed¶
Number of processed frames.
- Returns:
Number of processed frames.
- Return type:
int
- class pymolar.NdxFile(fname)¶
Bases:
objectReader for GROMACS NDX index files.
Example
ndx = pymolar.NdxFile("index.ndx") sel = ndx.get_group_as_sel("Protein", system)
- class pymolar.Particle¶
Bases:
objectView over one atom and its coordinates inside a system/selection.
Exposes both coordinate fields and atom descriptors as mutable properties.
Example
p = sys[0] print(p.name, p.resid, p.pos) # "CA", 1, [x, y, z] p.x = 1.5 # move to 1.5 nm on x axis p.mass = 12.0
- atomic_number¶
Atomic number.
- Returns:
Atomic number.
- Return type:
int
- bfactor¶
Temperature factor (B-factor).
- Returns:
B-factor value.
- Return type:
float
- chain¶
Chain identifier.
- Returns:
Chain identifier.
- Return type:
str
- charge¶
Atom charge.
- Returns:
Atom charge.
- Return type:
float
- id¶
Global atom index within the parent system/selection (read-only).
- Returns:
Global atom index.
- Return type:
int
- mass¶
Atomic mass.
- Returns:
Atomic mass.
- Return type:
float
- name¶
Atom name.
- Returns:
Atom name.
- Return type:
str
- occupancy¶
Occupancy value.
- Returns:
Occupancy value.
- Return type:
float
- pos¶
Return atom position as a length-3 NumPy array view.
- Returns:
Position vector
[x, y, z].- Return type:
numpy.ndarray
- resid¶
Residue identifier.
- Returns:
Residue identifier.
- Return type:
int
- resindex¶
Residue index.
- Returns:
Residue index.
- Return type:
int
- resname¶
Residue name.
- Returns:
Residue name.
- Return type:
str
- type_id¶
Force-field atom type identifier.
- Returns:
Force-field type id.
- Return type:
int
- type_name¶
Force-field atom type name.
- Returns:
Force-field type name.
- Return type:
str
- x¶
X coordinate.
- Returns:
X coordinate.
- Return type:
float
- y¶
Y coordinate.
- Returns:
Y coordinate.
- Return type:
float
- z¶
Z coordinate.
- Returns:
Z coordinate.
- Return type:
float
- class pymolar.PeriodicBox(*py_args)¶
Bases:
objectPeriodic simulation box geometry and minimum-image/PBC utilities.
Construct either from a 3x3 box matrix or from vectors and angles.
Example
import pymolar as molar import numpy as np box = molar.PeriodicBox(np.eye(3, dtype=np.float32) * 10.0) d = box.distance( np.array([0.0, 0.0, 0.0], dtype=np.float32), np.array([9.0, 0.0, 0.0], dtype=np.float32), [True, True, True], )
- closest_image(point, target, dims=None)¶
Return periodic image of
pointclosest totarget.- Parameters:
point – Point to image (length 3).
target – Reference target point (length 3).
dims – Periodic dimensions as
[x, y, z]booleans.
- Returns:
Imaged point closest to
target.- Return type:
numpy.ndarray
- Raises:
ValueError – If vector conversion fails.
- distance(p1, p2, dims)¶
Compute distance with optional periodic dimensions.
- Parameters:
p1 – First point (length 3).
p2 – Second point (length 3).
dims – Periodic dimensions as
[x, y, z]booleans.
- Returns:
Minimum-image distance in nm.
- Return type:
float
Example
d = box.distance([0.0, 0.0, 0.0], [1.0, 1.0, 1.0], [True, True, True])
- distance_squared(p1, p2, dims)¶
Compute squared distance with optional periodic dimensions.
- Parameters:
p1 – First point (length 3).
p2 – Second point (length 3).
dims – Periodic dimensions as
[x, y, z]booleans.
- Returns:
Squared minimum-image distance.
- Return type:
float
- Raises:
ValueError – If vector conversion fails.
- get_box_extents()¶
Return box extents in box basis.
- Returns:
Length-3 vector of extents.
- Return type:
numpy.ndarray
- get_lab_extents()¶
Return box extents in lab basis.
- Returns:
Length-3 vector of extents in lab basis.
- Return type:
numpy.ndarray
- get_matrix()¶
Return full box matrix.
- Returns:
3x3box matrix.- Return type:
numpy.ndarray
- is_triclinic()¶
Check whether the box is triclinic.
- Returns:
Truefor triclinic geometry, elseFalse.- Return type:
bool
- shortest_vector(arr, dims=None)¶
Compute minimum-image displacement vector.
- Parameters:
arr – Input displacement vector (length 3).
dims – Periodic dimensions as
[x, y, z]booleans.
- Returns:
Minimum-image displacement vector.
- Return type:
numpy.ndarray
- Raises:
ValueError – If vector conversion fails.
Example
v = box.shortest_vector([1.0, 0.0, 0.0])
- to_box_coords(point)¶
Convert Cartesian coordinates to box coordinates.
- Parameters:
point – Cartesian coordinate vector (length 3).
- Returns:
Coordinate in box basis.
- Return type:
numpy.ndarray
- Raises:
ValueError – If vector conversion fails.
- to_lab_coords(point)¶
Convert box coordinates to Cartesian coordinates.
- Parameters:
point – Coordinate in box basis (length 3).
- Returns:
Coordinate in Cartesian/lab basis.
- Return type:
numpy.ndarray
- Raises:
ValueError – If vector conversion fails.
- to_vectors_angles()¶
Return vector-length and angle representation.
- Returns:
Tuple
(vectors, angles)as NumPy arrays of length 3.- Return type:
tuple[numpy.ndarray, numpy.ndarray]
- wrap_point(p)¶
Wrap point into the primary unit cell.
- Parameters:
p – Input point (length 3).
- Returns:
Wrapped point.
- Return type:
numpy.ndarray
- Raises:
ValueError – If vector conversion fails.
Example
wrapped = box.wrap_point([5.0, 5.0, 5.0]) # back into unit cell
- class pymolar.Sasa¶
Bases:
objectSolvent-accessible surface area and volume measurements for a selection.
Both a result holder and a persistent calculator: call
update()on subsequent frames to reuse the power diagram without full reconstruction.Example
sasa = sel.sasa() print(sasa.total_area) # total SASA in nm² print(sasa.areas[:5]) # per-atom areas
- areas¶
Per-atom solvent-accessible areas (nm²).
- total_area¶
Total solvent-accessible area (nm²).
- total_volume¶
Total solvent-excluded volume (only meaningful when built with
sasa_vol).
- volumes¶
Per-atom solvent-excluded volumes (only meaningful when built with
sasa_vol).
- class pymolar.Sel¶
Bases:
objectAtom selection view with analysis and editing utilities.
Provides selection algebra, coordinate editing, and common analysis metrics.
Example
sel = system("name CA") com = sel.com() rg = sel.gyration()
- apply_transform(tr)¶
Apply rigid transform to selected coordinates.
- Parameters:
tr – Transform to apply.
- Returns:
None.- Return type:
None
Example
import pymolar tr = pymolar.fit_transform(mobile, ref) mobile.apply_transform(tr)
- box¶
Periodic box (proxied to backing state).
- Returns:
Periodic box.
- Return type:
- cog(dims=None)¶
Center of geometry, optionally using periodic dimensions.
- Parameters:
dims – Periodic dimensions
[x, y, z]booleans.- Returns:
Center-of-geometry vector
[x, y, z]in nm.- Return type:
numpy.ndarray
Example
center = sel.cog() # [x, y, z] in nm center = sel.cog(dims=[True, True, True]) # with PBC
- com(dims=None)¶
Center of mass, optionally using periodic dimensions.
- Parameters:
dims – Periodic dimensions
[x, y, z]booleans.- Returns:
Center-of-mass vector
[x, y, z]in nm.- Return type:
numpy.ndarray
Example
center = sel.com() # [x, y, z] in nm center = sel.com(dims=[True, True, True]) # with PBC
- coords¶
Coordinates of selected atoms as an array of shape
[3, n_atoms].- Returns:
Coordinate array.
- Return type:
numpy.ndarray
- gyration()¶
Radius of gyration (non-periodic).
- Returns:
Radius of gyration in nm.
- Return type:
float
Example
rg = sel.gyration() # radius of gyration in nm
- gyration_pbc()¶
Radius of gyration with periodic handling.
- Returns:
Radius of gyration.
- Return type:
float
- index¶
Global atom indices of this selection.
- Returns:
Index array.
- Return type:
numpy.ndarray
- inertia()¶
Principal moments and axes of inertia tensor.
- Returns:
Tuple
(moments, axes)where moments is length-3 and axes is 3×3.- Return type:
tuple[numpy.ndarray, numpy.ndarray]
Example
moments, axes = sel.inertia()
- inertia_pbc()¶
Principal moments and axes of inertia tensor with periodic handling.
- Returns:
Tuple
(moments, axes).- Return type:
tuple[numpy.ndarray, numpy.ndarray]
- iter_atoms()¶
Iterate over selected atoms.
- Returns:
Atom iterator.
- Return type:
- iter_pos()¶
Iterate over selected positions.
- Returns:
Position iterator.
- Return type:
- min_max()¶
Axis-aligned bounding-box min and max coordinates.
- Returns:
Tuple
(min_xyz, max_xyz)as NumPy arrays.- Return type:
tuple[numpy.ndarray, numpy.ndarray]
Example
lo, hi = sel.min_max() # two [x, y, z] arrays in nm
- principal_transform()¶
Principal-axes alignment transform (non-periodic).
- Returns:
Rigid transform.
- Return type:
IsometryTransform
- principal_transform_pbc()¶
Principal-axes alignment transform with periodic handling.
- Returns:
Rigid transform.
- Return type:
IsometryTransform
- replace_state_deep(st)¶
Replace state data in-place by swapping with a compatible state.
- Parameters:
st – Compatible state object.
- Returns:
None.- Return type:
None
- save(fname)¶
Save topology and selected state coordinates to file.
- Parameters:
fname – Output file path.
- Returns:
None.- Return type:
None
Example
sel.save("subset.pdb")
- set_box_from(sys)¶
Copy periodic box from a system.
- Parameters:
sys – Source system.
- Returns:
None.- Return type:
None
- set_same_bfactor(val)¶
Set B-factor for all selected atoms in-place.
- Parameters:
val – New B-factor value.
- set_same_chain(val)¶
Set chain ID for all selected atoms in-place.
- Parameters:
val – New chain character.
Example
sel.set_same_chain('A')
- set_same_mass(val)¶
Set atomic mass for all selected atoms in-place.
- Parameters:
val – New mass in Da.
- set_same_name(val)¶
Set atom name for all selected atoms in-place.
- Parameters:
val – New atom name.
- set_same_resid(val)¶
Set residue ID for all selected atoms in-place.
- Parameters:
val – New residue ID.
- set_same_resname(val)¶
Set residue name for all selected atoms in-place.
- Parameters:
val – New residue name.
Example
sel.set_same_resname('ALA')
- split_chain()¶
Split selection into sub-selections by chain ID.
- Returns:
List of per-chain selections.
- Return type:
list[Sel]
Example
chains = sel.split_chain() # list of Sel, one per chain
- split_molecule()¶
Split selection into molecular connected components.
- Returns:
List of per-molecule selections.
- Return type:
list[Sel]
Example
mols = sel.split_molecule() # list of Sel, one per molecule
- split_resindex()¶
Split selection into sub-selections by residue index.
- Returns:
List of per-residue selections.
- Return type:
list[Sel]
Example
residues = sel.split_resindex() # list of Sel, one per residue
- system¶
A
Systemsharing this selection’s topology and state.- Returns:
System view.
- Return type:
- time¶
Selection time value (proxied to backing state).
- Returns:
Time value.
- Return type:
float
- to_gromacs_ndx(name)¶
Export selection indices in GROMACS NDX group format.
- Parameters:
name – NDX group name.
- Returns:
NDX formatted group block.
- Return type:
str
Example
ndx_str = sel.to_gromacs_ndx("Protein")
- translate(arg)¶
Translate selected coordinates by vector.
- Parameters:
arg – Translation vector
[x, y, z]in nm.- Returns:
None.- Return type:
None
Example
import numpy as np sel.translate(np.array([0.5, 0.0, 0.0], dtype=np.float32)) # shift 0.5 nm along x
- whole_chains()¶
Expand selection to include all atoms in the same chains as any selected atom.
- Returns:
Expanded selection.
- Return type:
Example
ca = sys("name CA") chain = ca.whole_chains() # all atoms in those chains
- class pymolar.SelAtomIterator¶
Bases:
objectIterator over selected atoms.
- class pymolar.SelPosIterator¶
Bases:
objectIterator over selected atom positions.
- class pymolar.State¶
Bases:
objectCoordinate frame with simulation time and periodic box.
Stores coordinates, optional periodic box, and time for one frame.
Example
import pymolar fh = pymolar.FileHandler("traj.xtc", "r") st = fh.read_state() print(st.time) # simulation time in ps print(st.box) # PeriodicBox sys.state = st # hot-swap coordinates into system
- box¶
Periodic box of this state.
- Returns:
Periodic box.
- Return type:
- set_box_from(arg)¶
Copy periodic box from a
SystemorSel.- Parameters:
arg – Source object (
SystemorSel).- Returns:
None.- Return type:
None
- time¶
Simulation time value.
- Returns:
Frame time.
- Return type:
float
- class pymolar.SysAtomIterator¶
Bases:
objectIterator over system atoms.
- class pymolar.SysParticleIterator¶
Bases:
objectIterator over system particles (atom + position pairs).
- class pymolar.SysPosIterator¶
Bases:
objectIterator over system atom positions.
- class pymolar.System(*py_args)¶
Bases:
objectCoupled topology + state object used as the primary analysis source.
Example
import pymolar sys = pymolar.System("protein.pdb") sel = sys("protein and name CA") for particle in sys: print(particle.name, particle.pos) sys.save("out.pdb")
- append(*args)¶
Append atoms from a
Sel, selection expression, or(Atom, position)pair.- Parameters:
args – Either one selection-like argument or
(Atom, [x, y, z]).- Returns:
None.- Return type:
None
Example
import numpy as np sys.append(other_sel) # append a Sel sys.append(pymolar.Atom(), np.array([1.0, 2.0, 3.0], dtype=np.float32))
- box¶
Periodic box of the current frame; raises if box is absent.
- Returns:
Periodic box.
- Return type:
- iter_atoms()¶
Iterate over atoms.
- iter_pos()¶
Iterate over atom positions.
- remove(arg)¶
Remove atoms selected by argument (
Sel, query string, range, or indices).- Parameters:
arg – Selection expression,
Sel, range tuple, or index list.- Returns:
None.- Return type:
None
Example
sys.remove("resname HOH")
- replace_state_deep(st)¶
Swap internal state data with st when layouts are compatible.
- save(fname)¶
Save topology and current state to file.
- set_box_from(src)¶
Copy periodic box from another System or Sel.
- time¶
Simulation time of the current frame in picoseconds.
- Returns:
Frame time in ps.
- Return type:
float
- class pymolar.Topology¶
Bases:
objectMolecular topology container (atoms and connectivity).
Example
top = sys.topology print(len(top)) # number of atoms
Functions¶
- pymolar.distance_search(cutoff, data1, data2=None, dims=None)¶
Find atom pairs within cutoff distance between one or two selections.
- Parameters:
cutoff – Distance cutoff in nm, or
"vdw"for van der Waals radii sum.data1 – First selection.
data2 – Second selection (optional; self-search if omitted).
dims – Periodic dimensions
[x, y, z]booleans.
- Returns:
Tuple
(pairs, distances)—pairsis[N, 2]index array,distancesis length-N float array.- Return type:
tuple[numpy.ndarray, numpy.ndarray]
Example
import pymolar pairs, dists = pymolar.distance_search(0.35, sel1, sel2) pairs, dists = pymolar.distance_search(0.35, sel1, dims=[True, True, True])
- pymolar.fit_transform(sel1, sel2)¶
Compute rigid transform that best aligns
sel1ontosel2.- Parameters:
sel1 – Mobile selection.
sel2 – Reference selection.
- Returns:
Rigid transform.
- Return type:
IsometryTransform
Example
import pymolar tr = pymolar.fit_transform(mobile_sel, ref_sel) mobile_sel.apply_transform(tr)
- pymolar.fit_transform_matching(sel1, sel2)¶
Align selections by matching atom names before fitting.
- Parameters:
sel1 – Mobile selection.
sel2 – Reference selection.
- Returns:
Rigid transform computed on matched atoms.
- Return type:
IsometryTransform
Example
import pymolar tr = pymolar.fit_transform_matching(mobile_sel, ref_sel) mobile_sel.apply_transform(tr)
- pymolar.greeting()¶
Print library greeting and version banner.
- pymolar.rmsd_mw(sel1, sel2)¶
Compute mass-weighted RMSD between two selections.
- Parameters:
sel1 – First selection.
sel2 – Second selection.
- Returns:
Mass-weighted RMSD in nm.
- Return type:
float
Example
import pymolar rw = pymolar.rmsd_mw(sel1, sel2)
- pymolar.rmsd_py(sel1, sel2)¶
Compute RMSD between two selections.
- Parameters:
sel1 – First selection.
sel2 – Second selection.
- Returns:
RMSD in nm.
- Return type:
float
Example
import pymolar r = pymolar.rmsd(sel1, sel2)