metafalcon package

Submodules

metafalcon.atoms module

Constains some atomic constants needed by the dynamics program.

metafalcon.constants module

Contains several constants needed by the dynamics program.

metafalcon.constants.from_energyunit(unit)[source]
metafalcon.constants.to_energyunit(unit)[source]

metafalcon.cvcustom module

In order to use user-defined CVs, modify the following functions and include this file in your working directory.

metafalcon.cvcustom.customcv(symbols, coords, prm)[source]

Return value of user-defined CV for given symbols, coordinates and prm.

Parameters:
  • symbols (list of str) – atomic symbols of the molecule
  • coords (numpy array) – coordinates of the molecule (shape: (N, 3))
  • prm (dict) – parameters given in the CV section of meta-config.py
Returns:

custom – CV value

Return type:

number

Warning

Do not modify function name, arguments or returns.

metafalcon.cvcustom.dcustomcv(symbols, coords, prm)[source]

Return derivative of user_defined CV wrt coordinates for given symbols, coordinates and prm.

Parameters:
  • symbols (list of strings) – atomic symbols of the molecule
  • coords (numpy array) – coordinates of the molecule (shape: (N, 3))
  • prm (dict) – parameters given in the CV section of meta-config.py
Returns:

dcustom – derivative of CV wrt coordinates (shape: (Nat, 3))

Return type:

numpy array

Warning

Do not modify function name, arguments or returns.

metafalcon.initial module

Contains functions for import and generation of initial conditions for MD simulations.

metafalcon.initial.get_initial_random(symbols, coords, temp)[source]

Sample random initial conditions at a given temperature.

Parameters:
  • symbols (list of str) – list of atomic symbols of the molecule
  • coords (np.2darray) – atomic coordinates of the molecule (shape: (N, 3))
  • temp (float) – average temperature to be sampled
Returns:

  • symbols (list of str) – list of atomic symbols of the molecule
  • coords (np.2darray) – atomic coordinates of the molecule (shape: (N, 3))
  • velos (np.2darray) – atomic velocities of the molecule (shape: (N, 3))

metafalcon.initial.get_temperature(symbols, velos)[source]

Get the temperature corresponding to a set of atomic velocities.

Parameters:
  • symbols (list of str) – list of atomic symbols of the molecule
  • velos (np.2darray) – atomic velocities of the molecule (shape: (N, 3))
Returns:

temp – temperature

Return type:

float

metafalcon.initial.get_temperature_distribution(symbols, coords, temp, n=10000)[source]

Get the temperature distribution of a set of randomly sampled initial conditions.

Parameters:
  • symbols (list of str) – list of atomic symbols of the molecule
  • coords (np.2darray) – atomic coordinates of the molecule (shape: (N, 3))
  • temp (float) – average temperature to be sampled
  • n (int) – number of initial conditions to be sampled to determine the temperature distribution
Returns:

temps – list of temperatures that have been determined for a set of initial conditions

Return type:

list of float

metafalcon.initial.plot_temperature_distribution(temps)[source]

Plot a histogram of the sampled temperatures.

Parameters:temps (list of float) – list of temperatures that have been determined for a set of initial conditions
metafalcon.initial.get_initial_from_log(step, logdir='.')[source]

Get initial coordinates and velocities of a given time step from a previous calculation.

Parameters:
  • step (int) – time step from which from which the initial conditions should be taken
  • logdir (str) – path to the results directory that contains the respective logfile (default: ‘.’)
Returns:

  • symbols (list of str) – list of atomic symbols of the molecule
  • coords (np.2darray) – atomic coordinates of the molecule (shape: (N, 3))
  • velos (np.2darray) – atomic velocities of the molecule (shape: (N, 3))

metafalcon.initial.get_laststep_from_log(logdir='.')[source]

Determine the latest step number that has been performed in an existing dynamics simulation.

Parameters:logdir (str) – path to the results directory that contains the respective logfile (default: ‘.’)
Returns:step – latest step number
Return type:int
metafalcon.initial.read_initial()[source]

Read initial coordinates and velocities from files in xyz-format.

Returns:
  • symbols (list of str) – list of atomic symbols of the molecule
  • coords (np.2darray) – atomic coordinates of the molecule (shape: (N, 3))
  • velos (np.2darray) – atomic velocities of the molecule (shape: (N, 3))
metafalcon.initial.write_initial(symbols, coords, velos)[source]

Write initial coordinates and velocities in xyz-format for use in the dynamics program.

Parameters:
  • symbols (list of str) – list of atomic symbols of the molecule
  • coords (np.2darray) – atomic coordinates of the molecule (shape: (N, 3))
  • velos (np.2darray) – atomic velocities of the molecule (shape: (N, 3))
metafalcon.initial.generate_initial(config)[source]

Generate initial coordinates and velocities from the information given in the config dict.

Parameters:config (dict) – molecular dynamics configuration dictionary as read from ‘meta-config.json’

metafalcon.md module

Contains a class for molecular dynamics simulation runs.

class metafalcon.md.moleculardynamics[source]

Class to carry out molecular dynamics simulations using the velocity-Verlet algorithm.

Variables:
  • symbols (list of str) – atomic symbols of the molecule
  • coords (np.2darray) – cartesian coordinates of the molecule (shape: (N, 3))
  • velos (np.2darray) – cartesian velocities of the molecule (shape: (N, 3))
  • accel (np.2darray) – cartesian accelerations determined dy some quantum chemistry code (shape: (N, 3))
  • masses (np.array) – masses of the atoms in the molecule in atomic mass units (shape: (N))
  • energy (float) – electronic energy in the current time step
  • filendx (int) – index of the current dynamics-*.log file
  • restart (bool) – if True, the previous results are read in and extended (default: False)
  • start (int) – first dynamics step to start with when doing a restart
  • dometa (bool) – modify the forces in each the dynamics step by the metadynamics scheme (default: False)
  • updatemeta (bool) – if False, the existing metadynamics potential is not updated (default: True)
read_input()[source]

Read configuration parameters from meta-config.json.

get_cfg(config, key, default=None)[source]
set_coordinates(coords)[source]

Set coordinates attribute.

Parameters:coords (np.2darray) – cartesian coordinates of the molecule (shape: (N, 3))
set_velocities(velos)[source]

Set velocities attribute.

Parameters:velos (np.2darray) – cartesian velocities of the molecule (shape: (N, 3))
get_centerofmass(coords)[source]

Return center of mass of the given molecule.

Parameters:coords (np.2darray) – cartesian coordinates of the molecule (shape: (N, 3))
Returns:centerofmass – x-, y-, z-coordinates of the center of mass (shape: (3))
Return type:np.array
shift_to_centerofmass(coords)[source]

Translate molecule into its center of mass.

Parameters:coords (np.2darray) – cartesian coordinates of the molecule (shape: (N, 3))
Returns:newcoords – shifted cartesian coordinates (shape: (N, 3))
Return type:np.2darray
momentum(velos)[source]

Return the total momentum.

Parameters:velos (np.2darray) – cartesian velocities of the molecule (shape: (N, 3))
Returns:mom – x-, y-, z-components of the total momentum of the molecule (shape: (3))
Return type:np.array
eliminate_trans(velos, mom)[source]

Eliminate translation from velocities.

Parameters:
  • velos (np.2darray) – cartesian velocities of the molecule (shape: (N, 3))
  • mom (np.array) – x-, y-, z-components of the total momentum of the molecule (shape: (3))
Returns:

newvelos – new velocities with eliminated translation

Return type:

np.2darray

angular_momentum(coords, velos)[source]

Return the total angular momentum.

Parameters:
  • coords (np.2darray) – cartesian coordinates of the molecule (shape: (N, 3))
  • velos (np.2darray) – cartesian velocities of the molecule (shape: (N, 3))
Returns:

angmom – x-, y-, z-components of the total angular momentum of the molecule (shape: (3))

Return type:

np.array

moment_of_inertia(coords)[source]

Return the moment of inertia.

Parameters:coords (np.2darray) – cartesian coordinates of the molecule (shape: (N, 3))
Returns:inertia – x-, y-, z-components of the moment of inertia of the molecule (shape: (3))
Return type:np.array
get_angular_velocity(angmom, inertia)[source]

Return the angular velocity.

Parameters:
  • angmom (np.array) – x-, y-, z-components of the total angular momentum of the molecule (shape: (3))
  • inertia (np.array) – x-, y-, z-components of the moment of inertia of the molecule (shape: (3))
Returns:

angvel – x-, y-, z-components of the angular velocity of the molecule (shape: (3))

Return type:

np.array

eliminate_rot(velos, angvel, coords)[source]

Eliminate rotation from velocities.

Parameters:
  • velos (np.2darray) – cartesian velocities of the molecule (shape: (N, 3))
  • angvel (np.array) – x-, y-, z-components of the angular velocity of the molecule (shape: (3))
  • coords (np.2darray) – cartesian coordinates of the molecule (shape: (N, 3))
Returns:

newvelos – new velocities with eliminated rotation (shape: (N, 3))

Return type:

np.2darray

get_ekin(velos)[source]

Return the total kinetic energy.

Parameters:velos (np.2darray) – cartesian velocities of the molecule (shape: (N, 3))
Returns:ekin – kinetic energy
Return type:float
get_temperature(ekin)[source]

Return the temperature corresponding to the kinetic energy.

Parameters:ekin (float) – current kinetic energy
Returns:temp – current temperature
Return type:float
set_energy_and_accel()[source]

Use the QC-calculator to determine the current energy and accelerations.

initiate_trajectory()[source]

Initialize trajectory for a molecular dynamics simulation.

write_initial_step()[source]

Initialize the dynamics-*.log file.

write_step(step, ekin, epot, etot, conservation)[source]

Write results of current dynamics step to the dynamics-*.log file.

verlet()[source]

Perform molecular dynamics simulation using the velocity-Verlet algorithm.

metafalcon.meta module

Contains a class for the modification of MD forces by the metadynamics scheme.

class metafalcon.meta.metadynamics(symbols=None, file_extension=None, restart=False, update=True, use_threshold=False)[source]

Class for performing metadynamics simulations.

Parameters:
  • symbols (list of str, optional) – element symbols needed for calculation of coordination numbers (default: None)
  • file_extension (string, optional) – string to be appended to in- and output filenames to distinguish between multiple metadynamics objects in the same folder (default: None)
  • restart (bool, optional) – read gaussian centers from vgcentersfile to restart previous run (default: False)
  • update (bool, optional) – whether metadynamics potential should be updated during the simulation (default: True)
  • use_threshold (bool, optional) – use threshold instead of fixed tau_g (default: False)
Variables:
  • symbols (list of str) – element symbols of the molecule
  • file_extension (string, optional) – string to be appended to output filenames to distinguish between multiple metadynamics objects in the same folder
  • restart (bool) – read gaussian centers from vgcentersfile to restart previous run
  • update (bool) – whether metadynamics potential should be updated during the simulation
  • use_threshold (bool, optional) – use threshold instead of fixed tau_g (default: False)
  • vg_centers (list of float) – list of all gaussian centers added to vg
  • vg_steps (list of int) – step numbers, in which gaussians have been added to vg
  • cfg (dict) – configuration dictionary read from the input file
  • cvs (list) – list that containing the CV objects
read_input()[source]

Parse input from the input file

get_cfg(key, default=None)[source]
set_dftb_object(dftb)[source]

For Mulliken charges, set the DFTB.PES.PotentialEnergySurfaces object.

initialize_vg()[source]

Initialize vg and output files

set_s_and_ds()[source]

Calculate values of CV s and derivatives wrt coordinates ds_dr

get_height()[source]

Return the height of the gaussian in the current step.

get_vg()[source]

Return vg value in the current step

Returns:vg – value of the metadynamics potential in the current step
Return type:float
get_dvg_ds()[source]

Return dvg_ds in the current step

Returns:dvg_ds – derivative of the metadynamics potential wrt CV in the current step (shape: (Ncvs))
Return type:np.array
set_thresh_value(value)[source]

Set value for the threshold evaluation.

update_vg(step=0)[source]

Add new gaussian to vg every tau_g steps.

Parameters:step (int) – step number of the dynamics simulation (default: 0)
get_forces(coords, step)[source]

Check if vg should be updated and return metadynamics forces.

Parameters:
  • coords (np.array) – coordinates (shape: (N, 3))
  • step (int) – step number of the dynamics simulation
Returns:

forces – forces obtained from the current metadynamics potential (shape: (N, 3))

Return type:

np.array

get_coupled_energies_grads(coords, en1, en2, grad1, grad2, step=0)[source]

Couple the given energies and gradients by metadynamics potential.

Parameters:
  • coords (np.array) – coordinates (shape: (N, 3))
  • en1 (float) – lower state energy to be coupled
  • en2 (float) – upper state energy to be coupled
  • grad1 (np.array) – lower state energy gradient to be coupled (shape: (N, 3))
  • grad2 (np.array) – upper state energy gradient to be coupled (shape: (N, 3))
  • step (int, optional) – step number of the dynamics simulation needed for metadynamics wih tau_g (default: 0)
Returns:

  • new_en1 (float) – lower state energy after coupling
  • new_en2 (float) – upper state energy after coupling
  • new_grad1 (np.array) – lower state energy gradient after coupling (shape: (N, 3))
  • new_grad2 (np.array) – upper state energy gradient after coupling (shape: (N, 3))

metafalcon.wall module

Contains different analytic potential wall functions for application on a CV.

metafalcon.wall.polynomial_upper(x, s=1.0, offset=0.0, exp=2.0, **kwargs)[source]

upper quartic wall potential that is used for x values higher than offset

Parameters:
  • x (1darray or float) – x value(s) to be evaluated
  • s (float) – scaling factor (default: 1.0)
  • offset (float) – horizontal shift of the wall (default: 0.0)
  • exp (float) – exponent (default: 2.0)
Returns:

y – potential for the given x values

Return type:

1darray or float

metafalcon.wall.dpolynomial_upper(x, s=1.0, offset=0.0, exp=2.0, **kwargs)[source]

derivative of upper quartic wall potential that is used for x values higher than offset

Parameters:
  • x (1darray or float) – x value(s) to be evaluated
  • s (float) – scaling factor (default: 1.0)
  • offset (float) – horizontal shift of the wall (default: 0.0)
  • exp (float) – exponent (default: 2.0)
Returns:

y – potential for the given x values

Return type:

1darray or float

metafalcon.wall.polynomial_lower(x, s=1.0, offset=0.0, exp=2.0, **kwargs)[source]

lower quartic wall potential that is used for x values lower than offset

Parameters:
  • x (1darray or float) – x value(s) to be evaluated
  • s (float) – scaling factor (default: 1.0)
  • offset (float) – horizontal shift of the wall (default: 0.0)
  • exp (float) – exponent (default: 2.0)
Returns:

y – potential for the given x values

Return type:

1darray or float

metafalcon.wall.dpolynomial_lower(x, s=1.0, offset=0.0, exp=2.0, **kwargs)[source]

derivative of lower quartic wall potential that is used for x values higher than offset

Parameters:
  • x (1darray or float) – x value(s) to be evaluated
  • s (float) – scaling factor (default: 1.0)
  • offset (float) – horizontal shift of the wall (default: 0.0)
  • exp (float) – exponent (default: 2.0)
Returns:

y – potential for the given x values

Return type:

1darray or float

metafalcon.wall.log_upper(x, offset=0.0, c1=5.0, c2=5.0, **kwargs)[source]

upper logarithmic wall potential that is used for x values higher than offset

Parameters:
  • x (1darray or float) – x value(s) to be evaluated
  • offset (float) – horizontal shift of the wall (default: 0.0)
  • c1 (float) – weight factor in (kcal mol^-1)^-1 (default: 5.0)
  • c2 (float) – strongness factor in kcal mol^-1 (default: 5.0)
Returns:

y – potential for the given x values

Return type:

1darray or float

metafalcon.wall.dlog_upper(x, offset=0.0, c1=5.0, c2=5.0, **kwargs)[source]

derivative of upper logarithmic wall potential that is used for x values higher than offset

Parameters:
  • x (1darray or float) – x value(s) to be evaluated
  • offset (float) – horizontal shift of the wall (default: 0.0)
  • c1 (float) – weight factor in (kcal mol^-1)^-1 (default: 5.0)
  • c2 (float) – strongness factor in kcal mol^-1 (default: 5.0)
Returns:

y – potential for the given x values

Return type:

1darray or float

metafalcon.wall.log_lower(x, offset=0.0, c1=5.0, c2=5.0, **kwargs)[source]

upper logarithmic wall potential that is used for x values lower than offset

Parameters:
  • x (1darray or float) – x value(s) to be evaluated
  • offset (float) – horizontal shift of the wall (default: 0.0)
  • c1 (float) – weight factor in (kcal mol^-1)^-1 (default: 5.0)
  • c2 (float) – strongness factor in kcal mol^-1 (default: 5.0)
Returns:

y – potential for the given x values

Return type:

1darray or float

metafalcon.wall.dlog_lower(x, offset=0.0, c1=5.0, c2=5.0, **kwargs)[source]

derivative of lower logarithmic wall potential that is used for x values lower than offset

Parameters:
  • x (1darray or float) – x value(s) to be evaluated
  • offset (float) – horizontal shift of the wall (default: 0.0)
  • c1 (float) – weight factor in (kcal mol^-1)^-1 (default: 5.0)
  • c2 (float) – strongness factor in kcal mol^-1 (default: 5.0)
Returns:

y – potential for the given x values

Return type:

1darray or float

metafalcon.wall.plot_wall(wallfunc, xmin=0.0, xmax=1.0, **kwargs)[source]

Plot the wall in the given x-range in order check parameters.

Parameters:

metafalcon.xyz module

Contains useful functions for the work with xyz-files and cartesian coordinates.

metafalcon.xyz.read_xyz(filename)[source]

Read single structure from xyz-file and return symbols and coordinates

Parameters:filename (string) – path to xyz-file to be read
Returns:
  • symbols (list of strings) – element symbols
  • coords (numpy array (shape: (N, 3))) – coordinates in atomic units
metafalcon.xyz.read_trajectory(filename, first=None, last=None)[source]

Read trajectory from xyz-file and return symbols and coordinates

Parameters:
  • filename (string) – path to xyz-file to be read
  • first (int) – first frame to be read (default: 0)
  • last (int) – last frame to be read (default: -1)
Returns:

  • symbols (list of lists of strings) – element symbols
  • coords (list of numpy arrays (shape: (N, 3))) – coordinates in atomic units

metafalcon.xyz.write_trajectory(filename, symbols, coords, comments=[], fmt='%12.8f', mode='w')[source]

Write trajectory to xyz-file

Parameters:
  • filename (string) – path to xyz-file to be saved
  • symbols (list of strings) – element symbols
  • coords (list of numpy arrays (shape: (N, 3))) – coordinates in atomic units
  • comments (list of str) – comments to be inserted into the second line of every structure
  • fmt (str) – formatstring to be used for float numbers in the file
  • mode (str) – file opening mode (w for new file, a for append to existing; default: w
metafalcon.xyz.write_xyz(filename, symbols, coords, comment='', fmt='%12.8f', mode='w')[source]

Write single structure to xyz-file

Parameters:
  • filename (string) – path to xyz-file to be saved
  • symbols (list of strings) – element symbols
  • coords (numpy array (shape: (N, 3))) – coordinates in atomic units
  • comment (str) – comment to be inserted into the second line of the xyz-file
  • fmt (str) – formatstring to be used for float numbers in the file
  • mode (str) – file opening mode (w for new file, a for append to existing; default: w)
metafalcon.xyz.rot_x(coords, angle)[source]

Rotate molecule around x-axis

Parameters:
  • coords (numpy array (shape: (N, 3))) – coordinates of the molecule
  • angle (float) – angle in radians
Returns:

coords – rotated coordinates of the molecule

Return type:

numpy array (shape: (N, 3))

metafalcon.xyz.rot_y(coords, angle)[source]

Rotate molecule around y-axis

Parameters:
  • coords (numpy array (shape: (N, 3))) – coordinates of the molecule
  • angle (float) – angle in radians
Returns:

coords – rotated coordinates of the molecule

Return type:

numpy array (shape: (N, 3))

metafalcon.xyz.rot_z(coords, angle)[source]

Rotate molecule around z-axis

Parameters:
  • coords (numpy array (shape: (N, 3))) – coordinates of the molecule
  • angle (float) – angle in radians
Returns:

coords – rotated coordinates of the molecule

Return type:

numpy array (shape: (N, 3))

metafalcon.xyz.rot_rodrigues(coords, angle, n, center=array([0., 0., 0.]))[source]

Rotate molecule around arbitrary axis n

Parameters:
  • coords (numpy array (shape: (N, 3))) – coordinates of the molecule
  • angle (float) – angle in radians
  • n (numpy array (shape: 3)) – normed rotation vector
  • center (numpy array (shape: 3)) – rotation center (default: origin)
Returns:

coords – rotated coordinates of the molecule

Return type:

numpy array (shape: (N, 3))

metafalcon.xyz.trans_x(coords, t)[source]

Translate molecule in x-direction

Parameters:
  • coords (numpy array (shape: (N, 3))) – coordinates of the molecule
  • t (float) – translation
Returns:

coords – translated coordinates of the molecule

Return type:

numpy array (shape: (N, 3))

metafalcon.xyz.trans_y(coords, t)[source]

Translate molecule in y-direction

Parameters:
  • coords (numpy array (shape: (N, 3))) – coordinates of the molecule
  • t (float) – translation
Returns:

coords – translated coordinates of the molecule

Return type:

numpy array (shape: (N, 3))

metafalcon.xyz.trans_z(coords, t)[source]

Translate molecule in z-direction

Parameters:
  • coords (numpy array (shape: (N, 3))) – coordinates of the molecule
  • t (float) – translation
Returns:

coords – translated coordinates of the molecule

Return type:

numpy array (shape: (N, 3))

Module contents

Package for the input preparation, execution and analysis of molecular dynamics, metadynamics and multistate metadynamics simulations.