metafalcon.analysis package

Submodules

metafalcon.analysis.analyze module

Contains a class for the analysis of molecular dynamics trajectories.

class metafalcon.analysis.analyze.analyze_trajectory[source]

Bases: metafalcon.md.moleculardynamics

read_input()[source]

Read config dictionary from meta-config.json.

save_input()[source]

Save config dictionary to meta-config.json.

get_properties()[source]

Return all properties that have been specified to be analyzed.

add_property(ptype, **kwargs)[source]

Add a property to be analyzed along the calculated trajectory.

Parameters:
  • ptype (str) – property type to be analyzed (see below for all possibilities)
  • kwargs – see below
Keyword Arguments:
 
  • atoms (list) – atom indices that specify the atoms to be used for bond, angle or torsion analysis
  • atom (int) – atom index for the calculation of the coordination number
  • reference (list of int) – indices of reference atoms to be included (default: all atoms other than atom)
  • n (int) – first exponent
  • m (int) – second exponent
  • d ('auto' or float or list of float) – single distance or list of distances for all atoms in reference (default: ‘auto’)

Notes

The following properties are available for analysis:

Property types Explanation
pot potential energy
kin kinetic energy
temp temperature
cons energy conservation compared to last step (only NVE simulations)
bond bond length (requires atoms kwarg)
angle bond angle (requires atoms kwarg)
torsion torsion angle (requires atoms kwarg)
cn coordination number (requires atom, reference, m, n, and d kwargs)

Example

In order to analyze the bond length between the atoms 0 and 1, do:

>>> from metafalcon.analysis.analyze import analyze_trajectory
>>> anl = analyze_trajectory()
>>> anl.read_input()
>>> anl.add_property("bond", atoms=[0, 1])
>>> anl.analyze()
remove_property(ptype, **kwargs)[source]

Remove a property from the properties to be analyzed along the calculated trajectory. If the given parameters are not unique, all matching properties will be removed.

Parameters:
  • ptype (str) – property type to be removed (see below for all possibilities)
  • kwargs – see below
Keyword Arguments:
 
  • atoms (list) – atom indices that specify the atoms to be used for bond, angle or torsion analysis
  • atom (int) – atom index for the calculation of the coordination number
  • reference (list of int) – indices of reference atoms to be included (default: all atoms other than atom)
  • n (int) – first exponent
  • m (int) – second exponent
  • d ('auto' or float or list of float) – single distance or list of distances for all atoms in reference (default: ‘auto’)

Notes

The following properties are available for analysis:

Property types Explanation
pot potential energy
kin kinetic energy
temp temperature
cons energy conservation compared to last step (only NVE simulations)
bond bond length (requires atoms kwarg)
angle bond angle (requires atoms kwarg)
torsion torsion angle (requires atoms kwarg)
cn coordination number (requires atom, reference, m, n, and d kwargs)

Examples

Suppose you have added two bond lengths and a bond angle for analysis:

>>> from metafalcon.analysis.analyze import analyze_trajectory
>>> anl = analyze_trajectory()
>>> anl.read_input()
>>> anl.add_property("bond", atoms=[0, 1])
>>> anl.add_property("bond", atoms=[0, 2])
>>> anl.add_property("angle", atoms=[0, 1, 2])

You can remove a single property specifically

>>> anl.remove_property("bond", atoms=[0, 2])

or remove all bond lengths:

>>> anl.remove_property("bond")
reset_properties()[source]

Remove all properties that have been configured before.

set_stride(stride)[source]

Set the stride attribute to speed up the analysis.

Parameters:stride (int) – number of frames to be skipped after each analysis step (default: 1)
set_start(start)[source]

Set the first dynamics step to be analyzed.

Parameters:start (int) – first dynamics step to be analyzed (default: 0)
set_end(end)[source]

Set the last dynamics step to be analyzed.

Parameters:end (int) – last dynamics step to be analyzed (default: last)
set_nlog(nlog)[source]

Set the log-file to be analyzed.

Parameters:nlog (int) – log-file index to be analyzed
set_coordinates(coords)[source]

Set the coordinates attribute.

Parameters:coords (np.2darray) – molecular cartesian coordinates (shape: (N, 3))
log2xyz()[source]

Wrapper function for the coordinate conversion from the log-file to dynamics.xyz.

analyze()[source]

Wrapper function for the trajectory analysis of structure and energy properties.

metafalcon.analysis.reconstruct module

Reconstruct and plot the metadynamics potential.

class metafalcon.analysis.reconstruct.reconstruct_metadynamics(symbols=None, file_extension=None, restart=False, update=True, use_threshold=False)[source]

Bases: metafalcon.meta.metadynamics

reconstruct_vg(x0=None, ngauss=None, plot=True, animate=False, stride=1, unit='eV', fes=True)[source]

Reconstruct and plot potential from data in vg_centers.dat and meta-config.py. Grid and VG are saved as numpy arrays to vg_grid.npy and vg_values.npy

Parameters:
  • x0 (list of numpy arrays) – manually specified x-values for each collective variable (default: None)
  • ngauss (number) – first ngauss gaussian functions of the simulation will be evaluated (default: all)
  • plot (boolean) – plot single image of vg with the specified number of gaussians (default: True)
  • animate (boolean) – animate evolution of vg with number of gaussians (default: False)
  • stride (number) – number of gaussians to be added in each frame of the animation (default: 1)
  • fes (boolean) – plot the free energy surface instead of the metadynamics potential (default: True)
gauss(x, g, w)[source]
get_1dvg(x, heights, w, ngauss)[source]

Reconstruct vg on the 1D-grid x

Parameters:
  • x (numpy array) – x-values for the collective variable
  • h (number or list) – height of the added gaussians (single number for original, list for well-tempered)
  • w (number) – width of the added gaussians
  • ngauss (number) – first ngauss gaussian functions of the simulation will be evaluated
Returns:

vg – metadynamics potential on the grid x

Return type:

numpy array

get_ndvg(x, heights, w, ngauss)[source]

Reconstruct vg on the nD-grid x

Parameters:
  • x (list of numpy arrays) – x-values for each collective variable
  • h (number or list) – height of the added gaussians (single number for original, list for well-tempered)
  • w (list of numbers) – widths of the added gaussians (one per CV)
  • ngauss (number) – first ngauss gaussian functions of the simulation will be evaluated
Returns:

vg – metadynamics potential on the grid x

Return type:

numpy array

plot_1dvg(vg, x, cv, ngauss)[source]

Plot 1D metadynamics potential and save to file vg_ngauss-g.png

Parameters:
  • vg (1d numpy array) – metadynamics potential on the grid x
  • x (1d numpy array) – grid for the metadynamics potential
  • cv (cv object) – parameters of the collective variable read from input file
  • ngauss (number) – number of gaussians contained in vg
plot_2dvg(vg, x, cvs, ngauss)[source]

Plot 2D metadynamics potential and save to file vg_ngauss-g.png

Parameters:
  • vg (2d numpy array) – metadynamics potential on the grid x
  • x (2d numpy meshgrid) – grid for the metadynamics potential
  • cvs (list of cv objects) – list of parameters of the collective variables read from input file
  • ngauss (number) – number of gaussians contained in vg
animate_1dvg(vg, x, h, w, cv, ngauss, stride=1)[source]

Animate 1D metadynamics potential and save frames to directory animation

Parameters:
  • vg (1d numpy array) – metadynamics potential on the grid x
  • x (1d numpy array) – grid for the metadynamics potential
  • h (number or list) – height of the added gaussians (single number for original, list for well-tempered)
  • w (number) – width of the added gaussians
  • cv (cv object) – parameters of the collective variable read from input file
  • ngauss (number) – number of gaussians contained in vg
  • stride (number) – number of gaussians to be added in each frame of the animation (default: 1)
animate_2dvg(vg, x, h, w, cvs, ngauss, stride=1)[source]

Animate 2D metadynamics potential and save frames to directory animation

Parameters:
  • vg (2d numpy array) – metadynamics potential on the grid x
  • x (2d numpy array) – grid for the metadynamics potential
  • h (number or list) – height of the added gaussians (single number for original, list for well-tempered)
  • w (list of numbers) – widths of the added gaussians (1 per CV)
  • cvs (list of cv objects) – parameters of the collective variables read from input file
  • ngauss (number) – number of gaussians contained in vg
  • stride (number) – number of gaussians to be added in each frame of the animation (default: 1)

metafalcon.analysis.reset module

Reset dynamics output to a specified dynamics step in the dynamics-*.log-file.

metafalcon.analysis.reset.reset_log(step)[source]

Reset the dynamics-*.log file to the given step.

Parameters:step (int) – number of the last step to be kept in the log-file
metafalcon.analysis.reset.cap_file(maxlines, filename)[source]

Cap a file after the specified number of lines.

Parameters:
  • maxlines (int) – number of lines to be kept in the file
  • filename (str) – name of the file to be capped
metafalcon.analysis.reset.get_nat()[source]

Get the number of atoms from the dynamics-0.log file.

Returns:nat – number of atoms
Return type:int
metafalcon.analysis.reset.cap_gradfile(maxgrad, filename)[source]

Cap a gradient file after the specified number of gradient steps.

Parameters:
  • maxgrad (int) – number of gradient steps to be kept in the file
  • filename (str) – name of the file to be capped
metafalcon.analysis.reset.reset_vg(step, extension='')[source]

Reset the files vg_centers.dat and vg_steps.dat to the given step.

Parameters:
  • step (int) – number of time steps to be kept
  • extension (str) – filename extension of the metadynamics object
metafalcon.analysis.reset.add_restart()[source]

Add restart option to the meta-config.json file.

Module contents

Subpackage of metafalcon that provides several possibilities for the analysis of trajectories.