deepmd.infer#
Submodules#
Classes#
Functions#
| Python interface to calculate model deviation. |
| Factory function that forwards to DeepEval (for compatibility). |
Package Contents#
- class deepmd.infer.DeepEval(model_file: str, *args: Any, auto_batch_size: bool | int | deepmd.utils.batch_size.AutoBatchSize = True, neighbor_list: ase.neighborlist.NewPrimitiveNeighborList | None = None, **kwargs: 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_file
Path
The name of the frozen model file.
- *args
list
Positional arguments.
- auto_batch_sizebool or
int
orAutoBatchSize
, default:True
If True, automatic batch size will be used. If int, it will be used as the initial batch size.
- neighbor_list
ase.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.
- **kwargs
dict
Keyword arguments.
- model_file
- deep_eval#
- property output_def: deepmd.dpmodel.output_def.ModelOutputDef#
- Abstractmethod:
Returns the output variable definitions.
- _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: 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()
- 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.
- class deepmd.infer.DeepPot(model_file: str, *args: Any, auto_batch_size: bool | int | deepmd.utils.batch_size.AutoBatchSize = True, neighbor_list: ase.neighborlist.NewPrimitiveNeighborList | None = None, **kwargs: Any)[source]#
Bases:
deepmd.infer.deep_eval.DeepEval
Potential energy model.
- Parameters:
- model_file
Path
The name of the frozen model file.
- *args
list
Positional arguments.
- auto_batch_sizebool or
int
orAutoBatchSize
, default:True
If True, automatic batch size will be used. If int, it will be used as the initial batch size.
- neighbor_list
ase.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.
- **kwargs
dict
Keyword arguments.
- model_file
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: Literal[True], fparam: numpy.ndarray | None, aparam: numpy.ndarray | None, mixed_type: bool, **kwargs: Any) tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray, numpy.ndarray] [source]#
- eval(coords: numpy.ndarray, cells: numpy.ndarray | None, atom_types: list[int] | numpy.ndarray, atomic: Literal[False], fparam: numpy.ndarray | None, aparam: numpy.ndarray | None, mixed_type: bool, **kwargs: Any) tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray]
- eval(coords: numpy.ndarray, cells: numpy.ndarray | None, atom_types: list[int] | numpy.ndarray, atomic: bool, fparam: numpy.ndarray | None, aparam: numpy.ndarray | None, mixed_type: bool, **kwargs: Any) tuple[numpy.ndarray, Ellipsis]
Evaluate energy, force, and virial. If atomic is True, also return atomic energy and atomic virial.
- Parameters:
- coords
np.ndarray
The coordinates of the atoms, in shape (nframes, natoms, 3).
- cells
np.ndarray
The cell vectors of the system, in shape (nframes, 9). If the system is not periodic, set it to None.
- atom_types
list
[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.
- fparam
np.ndarray
,optional
The frame parameters, by default None.
- aparam
np.ndarray
,optional
The atomic parameters, by default None.
- mixed_typebool,
optional
Whether the atom_types is mixed type, by default False.
- **kwargs
dict
[str
,Any
] Keyword arguments.
- coords
- 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:
- coord
numpy.ndarray
, n_frames x n_atoms x 3 Coordinates of system to calculate
- box
numpy.ndarray
orNone
, n_frames x 3 x 3 Box to specify periodic boundary condition. If None, no pbc will be used
- atype
numpy.ndarray
, n_atoms x 1 Atom types
- models
list
of
DeepPot
models
Models used to evaluate deviation
- fname
str
orNone
File to dump results, default None
- frequency
int
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
- fparam
numpy.ndarray
frame specific parameters
- aparam
numpy.ndarray
atomic specific parameters
- real_data
dict
,optional
real data to calculate RMS real error
- atomicbool, default:
False
If True, calculate the force model deviation of each atom.
- relative
float
, 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_v
float
, 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.
- coord
- Returns:
- model_devi
numpy.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.
- model_devi
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 compatibility).
- Parameters:
- *args
positional arguments
- **kwargs
keyword arguments
- Returns:
DeepEval
potentials