deepmd.infer

Submodules

Package Contents

Classes

DeepEval

High-level Deep Evaluator interface.

DeepPot

Potential energy model.

Functions

calc_model_devi(coord, box, atype, models[, fname, ...])

Python interface to calculate model deviation.

DeepPotential(→ deep_eval.DeepEval)

Factory function that forwards to DeepEval (for compatbility).

class deepmd.infer.DeepEval(model_file: str, *args: List[Any], auto_batch_size: bool | int | deepmd.utils.batch_size.AutoBatchSize = True, neighbor_list: ase.neighborlist.NewPrimitiveNeighborList | None = None, **kwargs: Dict[str, Any])[source]

Bases: abc.ABC

High-level Deep Evaluator interface.

The specific DeepEval, such as DeepPot and DeepTensor, should inherit from this class. This class provides a high-level interface on the top of the low-level interface.

Parameters:
model_filePath

The name of the frozen model file.

*argslist

Positional arguments.

auto_batch_sizebool or int or AutoBatchSize, default: True

If True, automatic batch size will be used. If int, it will be used as the initial batch size.

neighbor_listase.neighborlist.NewPrimitiveNeighborList, optional

The ASE neighbor list class to produce the neighbor list. If None, the neighbor list will be built natively in the model.

**kwargsdict

Keyword arguments.

abstract property output_def: deepmd.dpmodel.output_def.ModelOutputDef

Returns the output variable definitions.

property has_efield: bool

Check if the model has efield.

property has_spin: bool

Check if the model has spin.

get_rcut() float[source]

Get the cutoff radius of this model.

get_ntypes() int[source]

Get the number of atom types of this model.

get_type_map() List[str][source]

Get the type map (element name of the atom types) of this model.

get_dim_fparam() int[source]

Get the number (dimension) of frame parameters of this DP.

get_dim_aparam() int[source]

Get the number (dimension) of atomic parameters of this DP.

_get_natoms_and_nframes(coords: numpy.ndarray, atom_types: numpy.ndarray, mixed_type: bool = False) Tuple[int, int][source]
_expande_atype(atype: numpy.ndarray, nframes: int, mixed_type: bool)[source]
eval_descriptor(coords: numpy.ndarray, cells: numpy.ndarray | None, atom_types: numpy.ndarray, fparam: numpy.ndarray | None = None, aparam: numpy.ndarray | None = None, mixed_type: bool = False, **kwargs: Dict[str, Any]) numpy.ndarray[source]

Evaluate descriptors by using this DP.

Parameters:
coords

The coordinates of atoms. The array should be of size nframes x natoms x 3

cells

The cell of the region. If None then non-PBC is assumed, otherwise using PBC. The array should be of size nframes x 9

atom_types

The atom types The list should contain natoms ints

fparam

The frame parameter. The array can be of size : - nframes x dim_fparam. - dim_fparam. Then all frames are assumed to be provided with the same fparam.

aparam

The atomic parameter The array can be of size : - nframes x natoms x dim_aparam. - natoms x dim_aparam. Then all frames are assumed to be provided with the same aparam. - dim_aparam. Then all frames and atoms are provided with the same aparam.

efield

The external field on atoms. The array should be of size nframes x natoms x 3

mixed_type

Whether to perform the mixed_type mode. If True, the input data has the mixed_type format (see doc/model/train_se_atten.md), in which frames in a system may have different natoms_vec(s), with the same nloc.

Returns:
descriptor

Descriptors.

eval_typeebd() numpy.ndarray[source]

Evaluate output of type embedding network by using this model.

Returns:
np.ndarray

The output of type embedding network. The shape is [ntypes, o_size], where ntypes is the number of types, and o_size is the number of nodes in the output layer.

Raises:
KeyError

If the model does not enable type embedding.

See also

deepmd.tf.utils.type_embed.TypeEmbedNet

The type embedding network.

Examples

Get the output of type embedding network of graph.pb:

>>> from deepmd.infer import DeepPotential
>>> dp = DeepPotential("graph.pb")
>>> dp.eval_typeebd()
_standard_input(coords, cells, atom_types, fparam, aparam, mixed_type)[source]
get_sel_type() List[int][source]

Get the selected atom types of this model.

Only atoms with selected atom types have atomic contribution to the result of the model. If returning an empty list, all atom types are selected.

_get_sel_natoms(atype) int[source]
get_ntypes_spin() int[source]

Get the number of spin atom types of this model. Only used in old implement.

class deepmd.infer.DeepPot(model_file: str, *args: List[Any], auto_batch_size: bool | int | deepmd.utils.batch_size.AutoBatchSize = True, neighbor_list: ase.neighborlist.NewPrimitiveNeighborList | None = None, **kwargs: Dict[str, Any])[source]

Bases: deepmd.infer.deep_eval.DeepEval

Potential energy model.

Parameters:
model_filePath

The name of the frozen model file.

*argslist

Positional arguments.

auto_batch_sizebool or int or AutoBatchSize, default: True

If True, automatic batch size will be used. If int, it will be used as the initial batch size.

neighbor_listase.neighborlist.NewPrimitiveNeighborList, optional

The ASE neighbor list class to produce the neighbor list. If None, the neighbor list will be built natively in the model.

**kwargsdict

Keyword arguments.

Examples

>>> from deepmd.infer import DeepPot
>>> import numpy as np
>>> dp = DeepPot("graph.pb")
>>> coord = np.array([[1, 0, 0], [0, 0, 1.5], [1, 0, 3]]).reshape([1, -1])
>>> cell = np.diag(10 * np.ones(3)).reshape([1, -1])
>>> atype = [1, 0, 1]
>>> e, f, v = dp.eval(coord, cell, atype)

where e, f and v are predicted energy, force and virial of the system, respectively.

property output_def: deepmd.dpmodel.output_def.ModelOutputDef

Get the output definition of this model.

property output_def_mag: deepmd.dpmodel.output_def.ModelOutputDef

Get the output definition of this model with magnetic parts.

eval(coords: numpy.ndarray, cells: numpy.ndarray | None, atom_types: List[int] | numpy.ndarray, atomic: bool = False, fparam: numpy.ndarray | None = None, aparam: numpy.ndarray | None = None, mixed_type: bool = False, **kwargs: Dict[str, Any]) Tuple[numpy.ndarray, Ellipsis][source]

Evaluate energy, force, and virial. If atomic is True, also return atomic energy and atomic virial.

Parameters:
coordsnp.ndarray

The coordinates of the atoms, in shape (nframes, natoms, 3).

cellsnp.ndarray

The cell vectors of the system, in shape (nframes, 9). If the system is not periodic, set it to None.

atom_typesList[int] or np.ndarray

The types of the atoms. If mixed_type is False, the shape is (natoms,); otherwise, the shape is (nframes, natoms).

atomicbool, optional

Whether to return atomic energy and atomic virial, by default False.

fparamnp.ndarray, optional

The frame parameters, by default None.

aparamnp.ndarray, optional

The atomic parameters, by default None.

mixed_typebool, optional

Whether the atom_types is mixed type, by default False.

**kwargsDict[str, Any]

Keyword arguments.

Returns:
energy

The energy of the system, in shape (nframes,).

force

The force of the system, in shape (nframes, natoms, 3).

virial

The virial of the system, in shape (nframes, 9).

atomic_energy

The atomic energy of the system, in shape (nframes, natoms). Only returned when atomic is True.

atomic_virial

The atomic virial of the system, in shape (nframes, natoms, 9). Only returned when atomic is True.

deepmd.infer.calc_model_devi(coord, box, atype, models, fname=None, frequency=1, mixed_type=False, fparam: numpy.ndarray | None = None, aparam: numpy.ndarray | None = None, real_data: dict | None = None, atomic: bool = False, relative: float | None = None, relative_v: float | None = None)[source]

Python interface to calculate model deviation.

Parameters:
coordnumpy.ndarray, n_frames x n_atoms x 3

Coordinates of system to calculate

boxnumpy.ndarray or None, n_frames x 3 x 3

Box to specify periodic boundary condition. If None, no pbc will be used

atypenumpy.ndarray, n_atoms x 1

Atom types

modelslist of DeepPot models

Models used to evaluate deviation

fnamestr or None

File to dump results, default None

frequencyint

Steps between frames (if the system is given by molecular dynamics engine), default 1

mixed_typebool

Whether the input atype is in mixed_type format or not

fparamnumpy.ndarray

frame specific parameters

aparamnumpy.ndarray

atomic specific parameters

real_datadict, optional

real data to calculate RMS real error

atomicbool, default: False

If True, calculate the force model deviation of each atom.

relativefloat, default: None

If given, calculate the relative model deviation of force. The value is the level parameter for computing the relative model deviation of the force.

relative_vfloat, default: None

If given, calculate the relative model deviation of virial. The value is the level parameter for computing the relative model deviation of the virial.

Returns:
model_devinumpy.ndarray, n_frames x 8

Model deviation results. The first column is index of steps, the other 7 columns are max_devi_v, min_devi_v, avg_devi_v, max_devi_f, min_devi_f, avg_devi_f, devi_e.

Examples

>>> from deepmd.tf.infer import calc_model_devi
>>> from deepmd.tf.infer import DeepPot as DP
>>> import numpy as np
>>> coord = np.array([[1, 0, 0], [0, 0, 1.5], [1, 0, 3]]).reshape([1, -1])
>>> cell = np.diag(10 * np.ones(3)).reshape([1, -1])
>>> atype = [1, 0, 1]
>>> graphs = [DP("graph.000.pb"), DP("graph.001.pb")]
>>> model_devi = calc_model_devi(coord, cell, atype, graphs)
deepmd.infer.DeepPotential(*args, **kwargs) deep_eval.DeepEval[source]

Factory function that forwards to DeepEval (for compatbility).

Parameters:
*args

positional arguments

**kwargs

keyword arguments

Returns:
DeepEval

potentials