Molssi Driver Interface Library
MDI Standard

Overview

The MDI Standard provides a straightforward, API-like method for enabling interoperability among computational molecular sciences codes. It uses a driver/engine model in which a driver code controls the high-level program flow of one or more engine codes. The driver exercises this control through the use of commands that are defined by the MDI Standard. Commands are available that correspond to a variety of tasks, such as "receive a new set of nuclear coordinates from me" (>COORDS), "start an MD simulation" (@INIT_MD) and "send me the forces on the nuclei" (<FORCES). The MDI standard defines the units, data types, and formatting of any data communicated between codes in response to a command.

Extending the MDI Standard

Although the MDI Standard has been designed to enable a wide range of potential applications, it is inevitable that some use cases will require functionality that is not yet supported by the standard. Typically, supporting such use cases will simply require addition of one or more new commands to the MDI Standard. Users who find themselves in need of a new MDI command are encouraged to create an issue describing their needs or submit a pull request suggesting a solution.

Units

All physical quantities communicated through MDI must be expressed in atomic units. Whenever it is necessary to perform a unit conversion on data communicated through MDI, the conversion factor should be obtained through the MDI_Conversion_Factor() function. Use of this function ensures self-consistency in unit conversions between codes, which is important for avoiding numerical instabilities.

Nodes

One of the powerful features of the MDI Standard is that it permits drivers to take advantage of existing implementations of time integrators and geometry optimizers in MDI engines. In particular, the @INIT_MD and @INIT_OPTG commands cause an engine to begin a molecular dynamics trajectory or a geometry optimization, respectively. Upon receiving one of these commands, an engine will perform the corresponding simulation without requiring further instruction from the driver, except that it will pause at certain "nodes" and listen for new commands.

The MDI Standard defines several nodes. Their names and when they occur are as follows:

@DEFAULT - Upon initially connecting to the driver.
@INIT_MC - Upon initializing a new Monte Carlo simulation.
@INIT_MD - Upon initializing a new molecular dynamics simulation.
@INIT_OPTG - Upon initializing a new geometry optimization.

Not all MDI engines are required to support each of these nodes, with the exception of the @DEFAULT node, which every MDI engine must support. Additional nodes may be defined by individual codes. Examples of possible code-specific nodes include:

@PRE-FORCES - After calculating all contributions to the atomic forces, except those associated with a constraint algorithm like SHAKE or RATTLE.
@FORCES - After calculating all contributions to the atomic forces.
@COORDS - After updating the atomic coordinates.

Several MDI commands enable a driver to control the program flow of an engine through different nodes. The <@ command instructs the engine to send the name of its node, while the @ command instructs the engine to proceed in its simulation until it reaches the next node. In addition, there are commands associated with each type of node (@PRE-FORCES, @FORCES, @COORDS) which instruct the engine to proceed in its simulation until it reaches the next node of that particular type.

A specific MD implementation might progress from a @PRE-FORCES node, to the @FORCES node, to the @COORDS node, and then repeat the cycle; however, this behavior is not guaranteed. Engines are permitted to pass through nodes in whatever order the implementation dictates. When writing drivers, it is a best practice to avoid assumptions about the ordering or frequency of nodes.

The MDI Library provides several functions that allow codes to define and communicate information about which nodes are supported:

Engines are not required to support every MDI command at every node. The MDI Library provides several functions that allow codes to define and communicate information about which commands are supported at specific nodes:

Constants

The following constants are defined by the MDI Standard and are accessible through the MDI Library:

  • MDI_MPI - Identifier for the MPI communication method
  • MDI_TCP - Identifier for the TCP communication method
  • MDI_LINK - Identifier for the LINK communication method
  • MDI_TEST - Identifier for the TEST communication method
  • MDI_DRIVER - Identifier for the DRIVER role
  • MDI_ENGINE - Identifier for the ENGINE role
  • MDI_INT - Data type identifier for integers
  • MDI_DOUBLE - Data type identifier for double precision floats
  • MDI_CHAR - Data type identifier for characters
  • MDI_INT_NUMPY - Data type identifier for Python NumPy integer arrays
  • MDI_DOUBLE_NUMPY - Data type identifier for Python NumPy double arrays
  • MDI_NAME_LENGTH - Maximum number of characters in the name of an MDI code (see Launching Codes with the MDI Library) or node (see Nodes)
  • MDI_COMMAND_LENGTH - Maximum number of characters in an MDI command (see Command List)

Command List

The following is a list of commands that are officially part of the MDI standard.

@

The engine proceeds to the next node (see Nodes). This command is typically not supported at the @DEFAULT node.

<@

The engine sends the driver a string that corresponds to the name of its current node (see Nodes).

Data Type: MDI_CHAR
Quantity: MDI_NAME_LENGTH

<CDENSITY

The engine sends the Cartesian coordinates of a set of grid points. This command is intended to be used in conjuction with the <NDENSITY and <DENSITY commands; these three commands enable a driver to acquire the electronic density distribution of an engine in a grid representation. See the <DENSITY command for more details.

Data Type: MDI_DOUBLE
Quantity: 3 * NDENSITY

>CELL

The driver sends a set of cell vectors to the engine, which resizes its simulation cell to the dimensions specified by the cell vectors.

Data Type: MDI_DOUBLE
Quantity: 9
Format: The first 3 values correspond to the x, y, and z values, respectively, of the first cell vector. The next 3 values correspond to the x, y, and z values, respectively, of the second cell vector. The next 3 values correspond to the x, y, and z values, respectively, of the third cell vector.

<CELL

The engine sends a set of cell vectors to the driver.

Data Type: MDI_DOUBLE
Quantity: 9
Format: The first 3 values correspond to the x, y, and z values, respectively, of the first cell vector. The next 3 values correspond to the x, y, and z values, respectively, of the second cell vector. The next 3 values correspond to the x, y, and z values, respectively, of the third cell vector.

>CELL_DISPL

The driver sends a displacement vector to the engine, which adjusts the origin of its simulation cell to the value of the displacement vector.

Data Type: MDI_DOUBLE
Quantity: 3
Format: The 3 values correspond to the x, y, and z values, respectively, of the simulation cell displacement vector.

<CELL_DISPL

The engine sends the displacement vector of the origin of its simulation cell to the driver.

Data Type: MDI_DOUBLE
Quantity: 3
Format: The 3 values correspond to the x, y, and z values, respectively, of the simulation cell displacement vector.

>CHARGES

The driver sends a set of atomic charges to the engine, which replaces its atomic charges with those sent by the driver.

Data Type: MDI_DOUBLE
Quantity: NATOMS
Format: Sequentially ascending order of atomic index

<CHARGES

The engine sends a set of atomic charges to the driver, in the same format as specified for the >CHARGES command.

Data Type: MDI_DOUBLE
Quantity: NATOMS
Format: Sequentially ascending order of atomic index

>CLATTICE

This command, along with the >NLATTICE and >LATTICE commands, allows the driver to assign a lattice of point charges to an engine, which incorporates the effects of these charges in all further calculations. After sending this command, the driver sends the coordinates of each of the point charges to the engine. Prior to sending this command, the driver must have set the number of point charges using the >NLATTICE command.

This command is primarily intended for use with gas-phase quantum mechanics codes. For an alternative command that is more appropriate for plane wave quantum mechanics codes, see the >POTENTIAL command.

Data Type: MDI_DOUBLE
Quantity: 3 * NLATTICE
Format: Sequentially ascending order of lattice charge index, with the coordinates for each individual lattice charge being provided in xyz order

@COORDS

The engine proceeds to the next @COORDS node (see Nodes). This command is not valid at the @DEFAULT node.

>COORDS

The driver sends a set of atomic coordinates to the engine, which replaces its atomic coordinates with those sent by the driver.

Data Type: MDI_DOUBLE
Quantity: 3 * NATOMS
Format: Sequentially ascending order of atomic index, with the coordinates for each individual atom being provided in xyz order

<COORDS

The engine sends a set of atomic coordinates to the driver.

Data Type: MDI_DOUBLE
Quantity: 3 * NATOMS
Format: Sequentially ascending order of atomic index, with the coordinates for each individual atom being provided in xyz order

>CPOTENTIAL

The driver sends the Cartesian coordinates of a set of grid points. This command is intended to be used in conjuction with the >NPOTENTIAL and >POTENTIAL commands; these three commands enable a driver to set an external potential that is incorporated into a subsequent scf_command command. See the >POTENTIAL command for more details.

Before sending this command, the driver must have first sent the number of grid points used to represent the potential via the >NPOTENTIAL command. It is also necessary that the driver send the values of the grid points via the >CPOTENTIAL command prior to any subsequent scf_command command.

Data Type: MDI_DOUBLE
Quantity: 3 * NPOTENTIAL

@DEFAULT

If not already at the @DEFAULT node, the engine exists whatever simulation (i.e. MD, OPTG, etc.) it is performing (possibly after completing an unfinished time step or geometry optimization step) and returns to the @DEFAULT node.

<DENSITY

The engine sends the value of its electronic density on a set of grid points. This command is intended to be used in conjuction with the <NDENSITY and <CDENSITY commands; these three commands enable a driver to acquire the electronic density distribution of an engine in a grid representation.

Data Type: MDI_DOUBLE
Quantity: NDENSITY

<DIMENSIONS

The engine sends basic information about the dimensionality of its system to the driver. For each of its three cell vectors (see the <CELL command) the engine sends an integer that indicates whether that dimension is represented as periodic, non-periodic, or not represented at all (in the case of 1d or 2d systems). The possible values for each cell vector are:

  • 0: Not represented
  • 1: Non-periodic
  • 2: Periodic
Data Type: MDI_INT
Quantity: 3
Format: Sequentially ascending order of cell vector (see the <CELL command)

>ELEC_MULT

The driver sends the electronic multiplicity of the system to the engine. This command is typically only appropriate for quantum mechanics engines.

Data Type: MDI_INT
Quantity: 1

<ELEC_MULT

The engine sends the electronic multiplicity of its system to the driver. This command is typically only appropriate for quantum mechanics engines.

Data Type: MDI_INT
Quantity: 1

<ELEMENTS

The engine sends the atomic number of each atom in its system to the driver.

Data Type: MDI_INT
Quantity: NATOMS
Format: Sequentially ascending order of atomic index

>ELEMENTS

The driver sends a set of atomic numbers to the engine, which replaces the atomic number of each atom in its system with the values sent by the driver.

Data Type: MDI_INT
Quantity: NATOMS
Format: Sequentially ascending order of atomic index

<ENERGY

If the engine is at the @DEFAULT node, it calculates and sends its total energy to the driver. If the engine has previously calculated the energy of the system, and no intervening commands from the driver could have changed the energy, the engine is permitted to send the previously calculated energy instead of recalculating it.

If the engine is not at the @DEFAULT node, it sends its most recently calculated total energy to the driver.

Data Type: MDI_DOUBLE
Quantity: 1

EXIT

The engine terminates and can no longer be sent commands.

@FORCES

The engine proceeds to the next @FORCES node (see Nodes). This command is not valid at the @DEFAULT node.

>FORCES

The driver sends a set of atomic forces to the engine, which replaces its internal forces with the forces sent by the driver.

Data Type: MDI_DOUBLE
Quantity: 3 * NATOMS
Format: Sequentially ascending order of atomic index, with the forces for each individual atom being provided in xyz order

>+FORCES

The driver sends a set of atomic forces to the engine, which adds the forces sent by the driver to its internal forces.

Data Type: MDI_DOUBLE
Quantity: 3 * NATOMS
Format: Sequentially ascending order of atomic index, with the forces for each individual atom being provided in xyz order

<FORCES

If the engine is at the @DEFAULT node, it calculates and sends its atomic forces to the driver. These forces include all force contributions, including the force contributions associated with any constraint algorithm (e.g. SHAKE, RATTLE, etc.). If the engine has previously calculated the atomic forces of the system, and no intervening commands from the driver could have changed the atomic forces, the engine is permitted to send the previously calculated atomic forces instead of recalculating them.

If the engine is not at the @DEFAULT node, it sends its most recently calculated atomic forces to the driver. Depending on the engine's current node, these forces may not include all contributions to the atomic forces. See the descriptions of the different nodes for more details.

Data Type: MDI_DOUBLE
Quantity: 3 * NATOMS
Format: Sequentially ascending order of atomic index, with the forces for each individual atom being provided in xyz order

@INIT_MC

The engine performs any initialization operations that are necessary before a Monte Carlo simulation can be performed, proceeding to the @INIT_MC node.

@INIT_MD

The engine performs any initialization operations that are necessary before a molecular dynamics simulation can be performed, proceeding to the @INIT_MD node.

Note: This command may change the engine's atomic coordinates under certain circumstances, such as if the SHAKE algorithm is used.

@INIT_OPTG

The engine performs any initialization operations that are necessary before a geometry optimization can be performed, proceeding to the @INIT_OPTG node.

Note: This command may change the engine's atomic coordinates under certain circumstances, such as if the SHAKE algorithm is used.

<KE

If the engine is at the @DEFAULT node, it calculates and sends its total kinetic energy to the driver. If the engine has previously calculated the energy of the system, and no intervening commands from the driver could have changed the energy, the engine is permitted to send the previously calculated energy instead of recalculating it.

If the engine is not at the @DEFAULT node, it sends its most recently calculated total kinetic energy to the driver.

Data Type: MDI_DOUBLE
Quantity: 1

<KE_ELEC

If the engine is at the @DEFAULT node, it calculates and sends the total kinetic energy of all electrons in its system to the driver. If the engine has previously calculated the energy of the system, and no intervening commands from the driver could have changed the energy, the engine is permitted to send the previously calculated energy instead of recalculating it.

If the engine is not at the @DEFAULT node, it sends its most recently calculated electronic kinetic energy to the driver.

Data Type: MDI_DOUBLE
Quantity: 1

<KE_NUC

If the engine is at the @DEFAULT node, it calculates and sends the total kinetic energy of all nuclei in its system to the driver. If the engine has previously calculated the energy of the system, and no intervening commands from the driver could have changed the energy, the engine is permitted to send the previously calculated energy instead of recalculating it.

If the engine is not at the @DEFAULT node, it sends its most recently calculated nuclear kinetic energy to the driver.

Data Type: MDI_DOUBLE
Quantity: 1

<LABELS

The engine sends a label for each atom in its system. "Labels" are intended primarily for the purpose of providing a human-readable identifier for each of the atoms, and do not have a standardized physical meaning. It is recommended that the labels correspond to the element of each atom (i.e., "H", "He", "Li", etc.), a name associated with atoms of a particular type (i.e., "Carboxyl_Hydrogen", "Methyl_Hydrogen"), or a similarly descriptive term. The atom labels may correspond to a number identifier (i.e., "1", "2", "3", etc.) in cases where more descriptive labels are not practical, but note that such labels must be represented using the MDI_CHAR data type, as indicated below. It is required that atoms having different physical properties (i.e., different force field terms in a molecular mechanics simulation or different nuclear charges in a quantum chemistry simulation) have different labels.

Data Type: MDI_CHAR
Quantity: MDI_LABEL_LENGTH * NATOMS
Format: An array of characters corresponding to the label of each atom in ascending order of atomic index, with each label consisting of MDI_LABEL_LENGTH characters and being padded with spaces (" ") where necessary

>LATTICE

This command, along with the >NLATTICE and >CLATTICE commands, allows the driver to assign a lattice of point charges to an engine, which incorporates the effects of these charges in all further calculations. After sending this command, the driver sends the charges of each of the point charges to the engine. Prior to sending this command, the driver must have set the number of point charges using the >NLATTICE command.

This command is primarily intended for use with gas-phase quantum mechanics codes. For an alternative command that is more appropriate for plane wave quantum mechanics codes, see the >POTENTIAL command.

Data Type: MDI_DOUBLE
Quantity: NLATTICE
Format: Sequentially ascending order of lattice charge index

<LATTICE_FORCES

If the engine is at the @DEFAULT node, it calculates and sends the forces on any lattice charges (which must have previously assigned with the >LATTICE command) to the driver. Prior to sending this command, the driver must have set the number, coordinates, and magnitudes of the lattice charges using the >NLATTICE, >CLATTICE, and >LATTICE commands. These forces must include only electrostatic interactions between the lattice charges and the atomic nuclei, and between the lattice charges and any electrons. They must not include electrostatic interactions between the lattice charges and other lattice charges. If the engine has previously calculated these forces, and no intervening commands from the driver could have changed the forces, the engine is permitted to send the previously calculated forces instead of recalculating them.

If the engine is not at the @DEFAULT node, it sends its most recently calculated lattice forces to the driver.

This command is primarily intended for use with gas-phase quantum mechanics codes.

Data Type: MDI_DOUBLE
Quantity: 3 * NLATTICE
Format: Sequentially ascending order of lattice charge index, with the forces for each individual lattice charge being provided in xyz order

<MASSES

The engine sends the driver the mass of each of the atoms.

Data Type: MDI_DOUBLE
Quantity: <NATOMS
Format: Sequentially ascending order of atomic index

>MONKHORST-PACK_NPOINTS

This command is typically expected for use with plane wave DFT engines. The driver sends the engine the number of k-points to generate on a Monkhorst-Pack grid. The engine then uses the k-points generated on this Monkhorst-Pack grid for all further simulations.

Data Type: MDI_INT
Quantity: 3
Format: The number of k-points to generate along each vector of the Brillouin zone, in ascending order of vector.

>MONKHORST-PACK_SHIFT

This command is typically expected for use with plane wave DFT engines. The driver sends the engine a set of values that indicate the extent to which a set of k-points on a Monkhorst-Pack grid should be displaced relative to the original (non-displaced) Monkhorst-Pack grid. The engine then uses the shifted k-points for all further simulations.

Data Type: MDI_DOUBLE
Quantity: 3
Format: The fraction of a grid step by which the k-points should be displaced, in ascending order of vector. A value of 0.0 indicates no displacement along the corresponding vector, while a value of 0.5 indicates a displacement of half a grid step in along the corresponding vector. Note that some engines can only support values of 0.0 or 0.5.

<NAME

The engine sends the driver a string that corresponds to the argument of -name in the MDI initialization options. This argument allows a driver to identify the purpose of connected engine codes within the simulation. For example, a particular QM/MM driver might require a connection with a single MM code and a single QM code, with the expected name of the MM code being "MM" and the expected name of the QM code being "QM". After initializing MDI and accepting communicators to the engines, the driver can use this command to identify which of the engines is the MM code and which is the QM code.

Data Type: MDI_CHAR
Quantity: MDI_NAME_LENGTH

<NATOMS

The engine sends the driver the number of atoms in the engine's system.

Data Type: MDI_INT
Quantity: 1

<NDENSITY

The engine sends the number of grid points it is using to represent its electronic density on a grid. This command is intended to be used in conjuction with the <CDENSITY and <DENSITY commands; these three commands enable a driver to acquire the electronic density distribution of an engine in a grid representation. See the <DENSITY command for more details.

Data Type: MDI_INT
Quantity: 1

>NLATTICE

This command, along with the >CLATTICE and >LATTICE commands, allows the driver to assign a lattice of point charges to an engine, which incorporates the effects of these charges in all further calculations. After sending this command, the driver sends the number of point charges to the engine. This command must be sent before either the >CLATTICE or >LATTICE commands can be sent.

This command is primarily intended for use with gas-phase quantum mechanics codes. For an alternative command that is more appropriate for plane wave quantum mechanics codes, see the >POTENTIAL command.

Data Type: MDI_INT
Quantity: 1

>NPOTENTIAL

The driver sends the number of grid points it is using to represent a potential on a grid. This command is intended to be used in conjuction with the >CPOTENTIAL and >POTENTIAL commands; these three commands enable a driver to set an external potential that is incorporated into a subsequent scf_command command. See the >POTENTIAL command for more details.

Data Type: MDI_INT
Quantity: 1

>POTENTIAL

The driver sends an set of values to the engine that correspond to a potential on a grid. If an scf_command command is later issued, this potential will be incorporated into the SCF calculation as an external potential.

Before sending this command, the driver must have first sent the number of grid points used to represent the potential via the >NPOTENTIAL command. It is also necessary that the driver send the Cartesian coordinates of the grid points via the >CPOTENTIAL command prior to any subsequent scf_command command.

Data Type: MDI_DOUBLE
Quantity: NPOTENTIAL

<PE

If the engine is at the @DEFAULT node, it calculates and sends its total potential energy to the driver. If the engine has previously calculated the energy of the system, and no intervening commands from the driver could have changed the energy, the engine is permitted to send the previously calculated energy instead of recalculating it.

If the engine is not at the @DEFAULT node, it sends its most recently calculated total potential energy to the driver.

Data Type: MDI_DOUBLE
Quantity: 1

<PE_ELEC

If the engine is at the @DEFAULT node, it calculates and sends its electronic potential energy to the driver. The electronic potential energy is defined as including all interactions between the electrons and any other particles or fields in the system. It also includes the interactions between the electrons and themselves. If the engine has previously calculated the energy of the system, and no intervening commands from the driver could have changed the energy, the engine is permitted to send the previously calculated energy instead of recalculating it.

If the engine is not at the @DEFAULT node, it sends its most recently calculated electronic potential energy to the driver.

Data Type: MDI_DOUBLE
Quantity: 1

<PE_NUC

If the engine is at the @DEFAULT node, it calculates and sends its nuclear potential energy to the driver. The nuclear potential energy is defined as including all interactions between the nuclei and any other particles or fields in the system, excluding any electrons. It also includes the interactions between the nuclei and themselves. If the engine has previously calculated the energy of the system, and no intervening commands from the driver could have changed the energy, the engine is permitted to send the previously calculated energy instead of recalculating it.

If the engine is not at the @DEFAULT node, it sends its most recently calculated electronic potential energy to the driver.

Data Type: MDI_DOUBLE
Quantity: 1

@PRE-FORCES

The engine proceeds to the next @PRE-FORCES node (see Nodes). This command is not valid at the @DEFAULT node.

<SPIN_POLARIZATION

The engine sends a value that indicates the manner in which it is currently simulating spin polarization. This command is typically intended for use with plane wave DFT engines.

Data Type: MDI_INT
Quantity: 1
Format: A value of 0 indicates that simulations will be performed in a non-spin-polarized manner. A value of 1 indicates that simulations will be performed in a spin-polarized manner, within the local spin density approximation (LSDA). A value of 2 indicates that simulations will be performed in a spin-polarized, noncollinear manner.

>SPIN_POLARIZATION

The driver sends a value that indicates the manner in which spin polarization should be simulated by the engine. This command is typically intended for use with plane wave DFT engines.

Data Type: MDI_INT
Quantity: 1
Format: A value of 0 indicates that simulations should be performed in a non-spin-polarized manner. A value of 1 indicates that simulations should be performed in a spin-polarized manner, within the local spin density approximation (LSDA). A value of 2 indicates that simulations should be performed in a spin-polarized, noncollinear manner.

<STRESS

The engine sends its virial stress tensor to the driver.

Data Type: MDI_DOUBLE
Quantity: 9
Format: The tensor components are sent in row-major order (xx, xy, xz, yx, yy, yz, zx, zy, zz).

>STRESS

The driver sends a virial stress tensor to the engine, which replaces its internal stress tensor with the stress tensor sent by the driver.

Data Type: MDI_DOUBLE
Quantity: 9
Format: The tensor components are sent in row-major order (xx, xy, xz, yx, yy, yz, zx, zy, zz).

>TOTCHARGE

The driver sends a value for the total charge of the system, including electron and nuclear charges, to the engine, which adjusts the number of electrons present in its system to the value required to reproduce the value sent by the driver. This command is typically only appropriate for quantum mechanics engines. Engines that support this command are not required to support non-integer charges; they are permitted to produce an error message if the value received deviates by more than 10^-12 from an integer, and to otherwise round the value received to the nearest integer.

Data Type: MDI_DOUBLE
Quantity: 1

<TOTCHARGE

The engine sends the total charge of its system, including electron and nuclear charges, to the driver.

Data Type: MDI_DOUBLE
Quantity: 1

>VELOCITES

The driver sends a set of atomic velocities to the driver, which replaces its atomic velocities with those provided by the driver.

Data Type: MDI_DOUBLE
Quantity: 3 * NATOMS
Format: Sequentially ascending order of atomic index, with the velocities for each individual atom being provided in xyz order

<VELOCITES

The engine sends the velocities of the atoms in its system to the driver.

Data Type: MDI_DOUBLE
Quantity: 3 * NATOMS
Format: Sequentially ascending order of atomic index, with the velocities for each individual atom being provided in xyz order

<VERSION

The engine sends the version number of the MDI Library to which it is linked to the driver.

Data Type: MDI_DOUBLE
Quantity: 1