deepmd package
Root of the deepmd package, exposes all public classes and submodules.
- class deepmd.DeepEval(model_file: Path, load_prefix: str = 'load', default_tf_graph: bool = False, auto_batch_size: Union[bool, int, deepmd.utils.batch_size.AutoBatchSize] = False)[source]
Bases:
object
Common methods for DeepPot, DeepWFC, DeepPolar, …
- Parameters
- model_file
Path
The name of the frozen model file.
- load_prefix: str
The prefix in the load computational graph
- default_tf_graphbool
If uses the default tf graph, otherwise build a new tf graph for evaluation
- auto_batch_sizebool or
int
orAutomaticBatchSize
, default:False
If True, automatic batch size will be used. If int, it will be used as the initial batch size.
- model_file
- Attributes
model_type
Get type of model.
model_version
Get version of model.
Methods
make_natoms_vec
(atom_types)Make the natom vector used by deepmd-kit.
reverse_map
(vec, imap)Reverse mapping of a vector according to the index map
sort_input
(coord, atom_type[, sel_atoms])Sort atoms in the system according their types.
- make_natoms_vec(atom_types: numpy.ndarray) numpy.ndarray [source]
Make the natom vector used by deepmd-kit.
- Parameters
- atom_types
The type of atoms
- Returns
natoms
The number of atoms. This tensor has the length of Ntypes + 2 natoms[0]: number of local atoms natoms[1]: total number of atoms held by this processor natoms[i]: 2 <= i < Ntypes+2, number of type i atoms
- static reverse_map(vec: numpy.ndarray, imap: List[int]) numpy.ndarray [source]
Reverse mapping of a vector according to the index map
- Parameters
- vec
Input vector. Be of shape [nframes, natoms, -1]
- imap
Index map. Be of shape [natoms]
- Returns
vec_out
Reverse mapped vector.
- static sort_input(coord: numpy.ndarray, atom_type: numpy.ndarray, sel_atoms: Optional[List[int]] = None)[source]
Sort atoms in the system according their types.
- Parameters
- coord
The coordinates of atoms. Should be of shape [nframes, natoms, 3]
- atom_type
The type of atoms Should be of shape [natoms]
- sel_atom
The selected atoms by type
- Returns
coord_out
The coordinates after sorting
atom_type_out
The atom types after sorting
idx_map
The index mapping from the input to the output. For example coord_out = coord[:,idx_map,:]
sel_atom_type
Only output if sel_atoms is not None The sorted selected atom types
sel_idx_map
Only output if sel_atoms is not None The index mapping from the selected atoms to sorted selected atoms.
- deepmd.DeepPotential(model_file: Union[str, pathlib.Path], load_prefix: str = 'load', default_tf_graph: bool = False) Union[deepmd.infer.deep_dipole.DeepDipole, deepmd.infer.deep_polar.DeepGlobalPolar, deepmd.infer.deep_polar.DeepPolar, deepmd.infer.deep_pot.DeepPot, deepmd.infer.deep_wfc.DeepWFC] [source]
Factory function that will inialize appropriate potential read from model_file.
- Parameters
- model_file: str
The name of the frozen model file.
- load_prefix: str
The prefix in the load computational graph
- default_tf_graphbool
If uses the default tf graph, otherwise build a new tf graph for evaluation
- Returns
Union
[DeepDipole
,DeepGlobalPolar
,DeepPolar
,DeepPot
,DeepWFC
]one of the available potentials
- Raises
RuntimeError
if model file does not correspond to any implementd potential
- class deepmd.DipoleChargeModifier(model_name: str, model_charge_map: List[float], sys_charge_map: List[float], ewald_h: float = 1, ewald_beta: float = 1)[source]
Bases:
deepmd.infer.deep_dipole.DeepDipole
- Parameters
- model_name
The model file for the DeepDipole model
- model_charge_map
Gives the amount of charge for the wfcc
- sys_charge_map
Gives the amount of charge for the real atoms
- ewald_h
Grid spacing of the reciprocal part of Ewald sum. Unit: A
- ewald_beta
Splitting parameter of the Ewald sum. Unit: A^{-1}
- Attributes
model_type
Get type of model.
model_version
Get version of model.
Methods
Build the computational graph for the force and virial inference.
eval
(coord, box, atype[, eval_fv])Evaluate the modification
eval_full
(coords, cells, atom_types[, ...])Evaluate the model with interface similar to the energy model.
get_dim_aparam
()Unsupported in this model.
get_dim_fparam
()Unsupported in this model.
get_ntypes
()Get the number of atom types of this model.
get_rcut
()Get the cut-off radius of this model.
get_sel_type
()Get the selected atom types of this model.
get_type_map
()Get the type map (element name of the atom types) of this model.
make_natoms_vec
(atom_types)Make the natom vector used by deepmd-kit.
modify_data
(data)Modify data.
reverse_map
(vec, imap)Reverse mapping of a vector according to the index map
sort_input
(coord, atom_type[, sel_atoms])Sort atoms in the system according their types.
- build_fv_graph() tensorflow.python.framework.ops.Tensor [source]
Build the computational graph for the force and virial inference.
- eval(coord: numpy.ndarray, box: numpy.ndarray, atype: numpy.ndarray, eval_fv: bool = True) Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray] [source]
Evaluate the modification
- Parameters
- coord
The coordinates of atoms
- box
The simulation region. PBC is assumed
- atype
The atom types
- eval_fv
Evaluate force and virial
- Returns
tot_e
The energy modification
tot_f
The force modification
tot_v
The virial modification
- modify_data(data: dict) None [source]
Modify data.
- Parameters
- data
Internal data of DeepmdData. Be a dict, has the following keys - coord coordinates - box simulation box - type atom types - find_energy tells if data has energy - find_force tells if data has force - find_virial tells if data has virial - energy energy - force force - virial virial
Subpackages
- deepmd.cluster package
- deepmd.descriptor package
- deepmd.entrypoints package
- Submodules
- deepmd.entrypoints.compress module
- deepmd.entrypoints.config module
- deepmd.entrypoints.convert module
- deepmd.entrypoints.doc module
- deepmd.entrypoints.freeze module
- deepmd.entrypoints.main module
- deepmd.entrypoints.test module
- deepmd.entrypoints.train module
- deepmd.entrypoints.transfer module
- deepmd.fit package
- deepmd.infer package
- deepmd.loggers package
- deepmd.loss package
- deepmd.model package
- deepmd.op package
- deepmd.utils package
- Submodules
- deepmd.utils.argcheck module
- deepmd.utils.batch_size module
- deepmd.utils.compat module
- deepmd.utils.convert module
- deepmd.utils.data module
- deepmd.utils.data_system module
- deepmd.utils.errors module
- deepmd.utils.graph module
- deepmd.utils.learning_rate module
- deepmd.utils.neighbor_stat module
- deepmd.utils.network module
- deepmd.utils.pair_tab module
- deepmd.utils.path module
- deepmd.utils.plugin module
- deepmd.utils.random module
- deepmd.utils.sess module
- deepmd.utils.tabulate module
- deepmd.utils.type_embed module
- deepmd.utils.weight_avg module
Submodules
deepmd.calculator module
ASE calculator interface module.
- class deepmd.calculator.DP(model: Union[str, pathlib.Path], label: str = 'DP', type_dict: Optional[Dict[str, int]] = None, **kwargs)[source]
Bases:
ase.calculators.calculator.Calculator
Implementation of ASE deepmd calculator.
Implemented propertie are energy, forces and stress
- Parameters
Examples
Compute potential energy
>>> from ase import Atoms >>> from deepmd.calculator import DP >>> water = Atoms('H2O', >>> positions=[(0.7601, 1.9270, 1), >>> (1.9575, 1, 1), >>> (1., 1., 1.)], >>> cell=[100, 100, 100], >>> calculator=DP(model="frozen_model.pb")) >>> print(water.get_potential_energy()) >>> print(water.get_forces())
Run BFGS structure optimization
>>> from ase.optimize import BFGS >>> dyn = BFGS(water) >>> dyn.run(fmax=1e-6) >>> print(water.get_positions())
- Attributes
- directory
- label
Methods
band_structure
()Create band-structure object for plotting.
calculate
([atoms, properties, system_changes])Run calculation with deepmd model.
calculate_numerical_forces
(atoms[, d])Calculate numerical forces using finite difference.
calculate_numerical_stress
(atoms[, d, voigt])Calculate numerical stress using finite difference.
calculate_properties
(atoms, properties)This method is experimental; currently for internal use.
check_state
(atoms[, tol])Check for any system changes since last calculation.
get_magnetic_moments
([atoms])Calculate magnetic moments projected onto atoms.
get_property
(name[, atoms, allow_calculation])Get the named property.
get_stresses
([atoms])the calculator should return intensive stresses, i.e., such that stresses.sum(axis=0) == stress
read
(label)Read atoms, parameters and calculated properties from output file.
reset
()Clear all information from old calculation.
set
(**kwargs)Set parameters like set(key1=value1, key2=value2, ...).
set_label
(label)Set label and convert label to directory and prefix.
calculation_required
export_properties
get_atoms
get_charges
get_default_parameters
get_dipole_moment
get_forces
get_magnetic_moment
get_potential_energies
get_potential_energy
get_stress
read_atoms
todict
- calculate(atoms: Optional[Atoms] = None, properties: List[str] = ['energy', 'forces', 'virial'], system_changes: List[str] = ['positions', 'numbers', 'cell', 'pbc', 'initial_charges', 'initial_magmoms'])[source]
Run calculation with deepmd model.
- Parameters
- atoms
Optional
[Atoms
],optional
atoms object to run the calculation on, by default None
- properties
List
[str
],optional
unused, only for function signature compatibility, by default [“energy”, “forces”, “stress”]
- system_changes
List
[str
],optional
unused, only for function signature compatibility, by default all_changes
- atoms
- implemented_properties: List[str] = ['energy', 'forces', 'virial', 'stress']
Properties calculator can handle (energy, forces, …)
- name = 'DP'
deepmd.common module
Collection of functions and classes used throughout the whole package.
- class deepmd.common.ClassArg[source]
Bases:
object
Class that take care of input json/yaml parsing.
The rules for parsing are defined by the add method, than parse is called to process the supplied dict
- Attributes
- arg_dict: Dict[str, Any]
dictionary containing parsing rules
- alias_map: Dict[str, Any]
dictionary with keyword aliases
Methods
add
(key, types_[, alias, default, must])Add key to be parsed.
get_dict
()Get dictionary built from rules defined by add method.
parse
(jdata)Parse input dictionary, use the rules defined by add method.
- add(key: str, types_: Union[type, List[type]], alias: Optional[Union[str, List[str]]] = None, default: Optional[Any] = None, must: bool = False) deepmd.common.ClassArg [source]
Add key to be parsed.
- Parameters
- Returns
ClassArg
instance with added key
- deepmd.common.add_data_requirement(key: str, ndof: int, atomic: bool = False, must: bool = False, high_prec: bool = False, type_sel: Optional[bool] = None, repeat: int = 1)[source]
Specify data requirements for training.
- Parameters
- key
str
type of data stored in corresponding *.npy file e.g. forces or energy
- ndof
int
number of the degrees of freedom, this is tied to atomic parameter e.g. forces have atomic=True and ndof=3
- atomicbool,
optional
specifies whwther the ndof keyworrd applies to per atom quantity or not, by default False
- mustbool,
optional
specifi if the *.npy data file must exist, by default False
- high_precbool,
optional
if tru load data to np.float64 else np.float32, by default False
- type_selbool,
optional
select only certain type of atoms, by default None
- repeat
int
,optional
if specify repaeat data repeat times, by default 1
- key
- deepmd.common.docstring_parameter(*sub: Tuple[str, ...])[source]
Add parameters to object docstring.
- Parameters
- sub: Tuple[str, …]
list of strings that will be inserted into prepared locations in docstring.
- deepmd.common.expand_sys_str(root_dir: Union[str, pathlib.Path]) List[str] [source]
Recursively iterate over directories taking those that contain type.raw file.
- deepmd.common.gelu(x: tensorflow.python.framework.ops.Tensor) tensorflow.python.framework.ops.Tensor [source]
Gaussian Error Linear Unit.
This is a smoother version of the RELU.
- Parameters
- x
tf.Tensor
float Tensor to perform activation
- x
- Returns
- x
with
the
GELU
activation
applied
- x
References
Original paper https://arxiv.org/abs/1606.08415
- deepmd.common.get_activation_func(activation_fn: _ACTIVATION) Callable[[tensorflow.python.framework.ops.Tensor], tensorflow.python.framework.ops.Tensor] [source]
Get activation function callable based on string name.
- Parameters
- activation_fn
_ACTIVATION
one of the defined activation functions
- activation_fn
- Returns
- Raises
RuntimeError
if unknown activation function is specified
- deepmd.common.get_np_precision(precision: _PRECISION) numpy.dtype [source]
Get numpy precision constant from string.
- Parameters
- precision
_PRECISION
string name of numpy constant or default
- precision
- Returns
np.dtype
numpy presicion constant
- Raises
RuntimeError
if string is invalid
- deepmd.common.get_precision(precision: _PRECISION) Any [source]
Convert str to TF DType constant.
- Parameters
- precision
_PRECISION
one of the allowed precisions
- precision
- Returns
tf.python.framework.dtypes.DType
appropriate TF constant
- Raises
RuntimeError
if supplied precision string does not have acorresponding TF constant
- deepmd.common.j_loader(filename: Union[str, pathlib.Path]) Dict[str, Any] [source]
Load yaml or json settings file.
- deepmd.common.j_must_have(jdata: Dict[str, _DICT_VAL], key: str, deprecated_key: List[str] = []) _DICT_VAL [source]
Assert that supplied dictionary conaines specified key.
- Returns
_DICT_VAL
value that was store unde supplied key
- Raises
RuntimeError
if the key is not present
- deepmd.common.make_default_mesh(test_box: numpy.ndarray, cell_size: float = 3.0) numpy.ndarray [source]
Get number of cells of size=`cell_size` fit into average box.
- Parameters
- test_box
np.ndarray
numpy array with cells of shape Nx9
- cell_size
float
,optional
length of one cell, by default 3.0
- test_box
- Returns
np.ndarray
mesh for supplied boxes, how many cells fit in each direction
- deepmd.common.select_idx_map(atom_types: numpy.ndarray, select_types: numpy.ndarray) numpy.ndarray [source]
Build map of indices for element supplied element types from all atoms list.
- Parameters
- atom_types
np.ndarray
array specifing type for each atoms as integer
- select_types
np.ndarray
types of atoms you want to find indices for
- atom_types
- Returns
np.ndarray
indices of types of atoms defined by select_types in atom_types array
Warning
select_types array will be sorted before finding indices in atom_types
deepmd.env module
Module that sets tensorflow working environment and exports inportant constants.
- deepmd.env.GLOBAL_ENER_FLOAT_PRECISION
alias of
numpy.float64
- deepmd.env.GLOBAL_NP_FLOAT_PRECISION
alias of
numpy.float64
- deepmd.env.global_cvt_2_ener_float(xx: tensorflow.python.framework.ops.Tensor) tensorflow.python.framework.ops.Tensor [source]
Cast tensor to globally set energy precision.