![]() |
ForceBalance API
1.3
Automated optimization of force fields and empirical potentials
|
Derived from Engine object for carrying out general purpose GROMACS calculations. More...
Public Member Functions | |
| def | __init__ (self, name="gmx", kwargs) |
| def | setopts (self, kwargs) |
| Called by init ; Set GROMACS-specific options. More... | |
| def | readsrc (self, kwargs) |
| Called by init ; read files from the source directory. More... | |
| def | prepare (self, pbc=False, kwargs) |
| Called by init ; prepare the temp directory and figure out the topology. More... | |
| def | get_charges (self) |
| def | links (self) |
| def | callgmx (self, command, stdin=None, print_to_screen=False, print_command=False, kwargs) |
| Call GROMACS; prepend the gmxpath to the call to the GROMACS program. More... | |
| def | warngmx (self, command, warnings=[], maxwarn=1, kwargs) |
| Call gromacs and allow for certain expected warnings. More... | |
| def | energy_termnames (self, edrfile=None) |
| Get a list of energy term names from the .edr file by parsing a system call to g_energy. More... | |
| def | optimize (self, shot, crit=1e-4, align=True, kwargs) |
| Optimize the geometry and align the optimized geometry to the starting geometry. More... | |
| def | evaluate_ (self, force=False, dipole=False, traj=None) |
| Utility function for computing energy, and (optionally) forces and dipoles using GROMACS. More... | |
| def | evaluate_snapshot (self, shot, force=False, dipole=False) |
| Evaluate variables (energies, force and/or dipole) using GROMACS for a single snapshot. More... | |
| def | evaluate_trajectory (self, force=False, dipole=False, traj=None) |
| Evaluate variables (energies, force and/or dipole) using GROMACS over a trajectory. More... | |
| def | make_gro_trajectory (self, fout=None) |
| Return the MD trajectory as a Molecule object. More... | |
| def | energy_one (self, shot) |
| Compute the energy using GROMACS for a snapshot. More... | |
| def | energy_force_one (self, shot) |
| Compute the energy and force using GROMACS for a single snapshot; interfaces with AbInitio target. More... | |
| def | energy (self, traj=None) |
| Compute the energy using GROMACS over a trajectory. More... | |
| def | energy_force (self, force=True, traj=None) |
| Compute the energy and force using GROMACS over a trajectory. More... | |
| def | energy_dipole (self, traj=None) |
| def | energy_rmsd (self, shot, optimize=True) |
| Calculate energy of the selected structure (optionally minimize and return the minimized energy and RMSD). More... | |
| def | interaction_energy (self, fraga, fragb) |
| Computes the interaction energy between two fragments over a trajectory. More... | |
| def | multipole_moments (self, shot=0, optimize=True, polarizability=False) |
| Return the multipole moments of the 1st snapshot in Debye and Buckingham units. More... | |
| def | normal_modes (self, shot=0, optimize=True) |
| def | generate_positions (self) |
| def | n_snaps (self, nsteps, step_interval, timestep) |
| def | scd_persnap (self, ndx, timestep, final_frame) |
| def | calc_scd (self, n_snap, timestep) |
| def | n_nonwater (self, structure_file) |
| def | molecular_dynamics (self, nsteps, timestep, temperature=None, pressure=None, nequil=0, nsave=0, minimize=True, threads=None, verbose=False, bilayer=False, kwargs) |
| Method for running a molecular dynamics simulation. More... | |
| def | md (self, nsteps=0, nequil=0, verbose=False, deffnm=None, kwargs) |
| Method for running a molecular dynamics simulation. More... | |
Public Attributes | |
| valkwd | |
| Valid GROMACS-specific keywords. More... | |
| gmxsuffix | |
| Disable some optimizations. More... | |
| gmx_eq_barostat | |
| Barostat keyword for equilibration. More... | |
| gmxpath | |
| The directory containing GROMACS executables (e.g. More... | |
| gmxversion | |
| Figure out the GROMACS version - it will determine how programs are called. More... | |
| top | |
| Attempt to determine file names of .gro, .top, and .mdp files. More... | |
| mdp | |
| mol | |
| gmx_defs | |
| pbc | |
| have_constraints | |
| Link files into the temp directory. More... | |
| double | |
| Write out the trajectory coordinates to a .gro file. More... | |
| AtomMask | |
| AtomLists | |
| mdtraj | |
| mdene | |
Derived from Engine object for carrying out general purpose GROMACS calculations.
| def src.gmxio.GMX.calc_scd | ( | self, | |
| n_snap, | |||
| timestep | |||
| ) |
| def src.gmxio.GMX.callgmx | ( | self, | |
| command, | |||
stdin = None, |
|||
print_to_screen = False, |
|||
print_command = False, |
|||
| kwargs | |||
| ) |
| def src.gmxio.GMX.energy | ( | self, | |
traj = None |
|||
| ) |
| def src.gmxio.GMX.energy_dipole | ( | self, | |
traj = None |
|||
| ) |
| def src.gmxio.GMX.energy_force | ( | self, | |
force = True, |
|||
traj = None |
|||
| ) |
| def src.gmxio.GMX.energy_force_one | ( | self, | |
| shot | |||
| ) |
| def src.gmxio.GMX.energy_one | ( | self, | |
| shot | |||
| ) |
| def src.gmxio.GMX.energy_rmsd | ( | self, | |
| shot, | |||
optimize = True |
|||
| ) |
| def src.gmxio.GMX.energy_termnames | ( | self, | |
edrfile = None |
|||
| ) |
| def src.gmxio.GMX.evaluate_ | ( | self, | |
force = False, |
|||
dipole = False, |
|||
traj = None |
|||
| ) |
Utility function for computing energy, and (optionally) forces and dipoles using GROMACS.
Inputs: force: Switch for calculating the force. dipole: Switch for calculating the dipole. traj: Trajectory file name. If present, will loop over these snapshots. Otherwise will do a single point evaluation at the current geometry.
Outputs: Result: Dictionary containing energies, forces and/or dipoles.
Definition at line 935 of file gmxio.py.
| def src.gmxio.GMX.evaluate_snapshot | ( | self, | |
| shot, | |||
force = False, |
|||
dipole = False |
|||
| ) |
| def src.gmxio.GMX.evaluate_trajectory | ( | self, | |
force = False, |
|||
dipole = False, |
|||
traj = None |
|||
| ) |
| def src.gmxio.GMX.generate_positions | ( | self | ) |
| def src.gmxio.GMX.get_charges | ( | self | ) |
| def src.gmxio.GMX.interaction_energy | ( | self, | |
| fraga, | |||
| fragb | |||
| ) |
| def src.gmxio.GMX.links | ( | self | ) |
| def src.gmxio.GMX.make_gro_trajectory | ( | self, | |
fout = None |
|||
| ) |
| def src.gmxio.GMX.md | ( | self, | |
nsteps = 0, |
|||
nequil = 0, |
|||
verbose = False, |
|||
deffnm = None, |
|||
| kwargs | |||
| ) |
Method for running a molecular dynamics simulation.
A little different than molecular_dynamics (for thermo.py)
Required arguments: nsteps = (int) Number of total time steps nequil = (int) Number of additional time steps at the beginning for equilibration verbose = (bool) Be loud and noisy deffnm = (string) default names for simulation output files The simulation data is written to the working directory.
Definition at line 1406 of file gmxio.py.
| def src.gmxio.GMX.molecular_dynamics | ( | self, | |
| nsteps, | |||
| timestep, | |||
temperature = None, |
|||
pressure = None, |
|||
nequil = 0, |
|||
nsave = 0, |
|||
minimize = True, |
|||
threads = None, |
|||
verbose = False, |
|||
bilayer = False, |
|||
| kwargs | |||
| ) |
Method for running a molecular dynamics simulation.
Required arguments: nsteps = (int) Number of total time steps timestep = (float) Time step in FEMTOSECONDS temperature = (float) Temperature control (Kelvin) pressure = (float) Pressure control (atmospheres) nequil = (int) Number of additional time steps at the beginning for equilibration nsave = (int) Step interval for saving data minimize = (bool) Perform an energy minimization prior to dynamics threads = (int) Number of MPI-threads
Returns simulation data: Rhos = (array) Density in kilogram m^-3 Potentials = (array) Potential energies Kinetics = (array) Kinetic energies Volumes = (array) Box volumes Dips = (3xN array) Dipole moments EComps = (dict) Energy components Als = (array) Average area per lipid in nm^2 Scds = (Nx28 array) Deuterium order parameter
Definition at line 1241 of file gmxio.py.
| def src.gmxio.GMX.multipole_moments | ( | self, | |
shot = 0, |
|||
optimize = True, |
|||
polarizability = False |
|||
| ) |
| def src.gmxio.GMX.n_snaps | ( | self, | |
| nsteps, | |||
| step_interval, | |||
| timestep | |||
| ) |
| def src.gmxio.GMX.normal_modes | ( | self, | |
shot = 0, |
|||
optimize = True |
|||
| ) |
| def src.gmxio.GMX.optimize | ( | self, | |
| shot, | |||
crit = 1e-4, |
|||
align = True, |
|||
| kwargs | |||
| ) |
| def src.gmxio.GMX.prepare | ( | self, | |
pbc = False, |
|||
| kwargs | |||
| ) |
| def src.gmxio.GMX.readsrc | ( | self, | |
| kwargs | |||
| ) |
| def src.gmxio.GMX.scd_persnap | ( | self, | |
| ndx, | |||
| timestep, | |||
| final_frame | |||
| ) |
| def src.gmxio.GMX.setopts | ( | self, | |
| kwargs | |||
| ) |
| def src.gmxio.GMX.warngmx | ( | self, | |
| command, | |||
warnings = [], |
|||
maxwarn = 1, |
|||
| kwargs | |||
| ) |
| src.gmxio.GMX.double |
Write out the trajectory coordinates to a .gro file.
At this point, we could have gotten a .mdp file from the target folder or as part of the force field. If it still missing, then we may write a default. Call grompp followed by gmxdump to read the trajectory
| src.gmxio.GMX.gmx_eq_barostat |
| src.gmxio.GMX.gmxpath |
| src.gmxio.GMX.gmxsuffix |
| src.gmxio.GMX.gmxversion |
Figure out the GROMACS version - it will determine how programs are called.
Always, always remove backup files.
Create symbolic links (mainly for the case of .top and .mdp files which don't exist at object creation) Call a GROMACS program as you would from the command line.
| src.gmxio.GMX.have_constraints |
| src.gmxio.GMX.top |
1.8.13