matscipy.calculators.eam package

Implements the Embedded Atom Method

Modules

matscipy.calculators.eam.calculator

EAM calculator

matscipy.calculators.eam.io

Read and write tabulated EAM potentials

matscipy.calculators.eam.average_atom

matscipy.calculators.eam.calculator module

EAM calculator

class matscipy.calculators.eam.calculator.EAM(fn=None, atomic_numbers=None, F=None, f=None, rep=None, cutoff=None, kind='eam/alloy')

Bases: MatscipyCalculator

implemented_properties: List[str] = ['energy', 'free_energy', 'stress', 'forces', 'hessian']

Properties calculator can handle (energy, forces, …)

default_parameters: Dict[str, Any] = {}

Default parameters

name = 'EAM'
__init__(fn=None, atomic_numbers=None, F=None, f=None, rep=None, cutoff=None, kind='eam/alloy')

Basic calculator implementation.

restart: str

Prefix for restart file. May contain a directory. Default is None: don’t restart.

ignore_bad_restart_file: bool

Deprecated, please do not use. Passing more than one positional argument to Calculator() is deprecated and will stop working in the future. Ignore broken or missing restart file. By default, it is an error if the restart file is missing or broken.

directory: str or PurePath

Working directory in which to read and write files and perform calculations.

label: str

Name used for all files. Not supported by all calculators. May contain a directory, but please use the directory parameter for that instead.

atoms: Atoms object

Optional Atoms object to which the calculator will be attached. When restarting, atoms will get its positions and unit-cell updated from file.

energy_virial_and_forces(atomic_numbers_i, i_n, j_n, dr_nc, abs_dr_n)

Compute the potential energy, the virial and the forces.

Parameters:
  • atomic_numbers_i (array_like) – Atomic number for each atom in the system

  • i_n (array_like) – Neighbor pairs

  • j_n (array_like) – Neighbor pairs

  • dr_nc (array_like) – Distance vectors between neighbors

  • abd_dr_n (array_like) – Length of distance vectors between neighbors

Returns:

  • epot (float) – Potential energy

  • virial_v (array) – Virial

  • forces_ic (array) – Forces acting on each atom

calculate(atoms, properties, system_changes)

Do the calculation.

properties: list of str

List of what needs to be calculated. Can be any combination of ‘energy’, ‘forces’, ‘stress’, ‘dipole’, ‘charges’, ‘magmom’ and ‘magmoms’.

system_changes: list of str

List of what has changed since last calculation. Can be any combination of these six: ‘positions’, ‘numbers’, ‘cell’, ‘pbc’, ‘initial_charges’ and ‘initial_magmoms’.

Subclasses need to implement this, but can ignore properties and system_changes if they want. Calculated properties should be inserted into results dictionary like shown in this dummy example:

self.results = {'energy': 0.0,
                'forces': np.zeros((len(atoms), 3)),
                'stress': np.zeros(6),
                'dipole': np.zeros(3),
                'charges': np.zeros(len(atoms)),
                'magmom': 0.0,
                'magmoms': np.zeros(len(atoms))}

The subclass implementation should first call this implementation to set the atoms attribute and create any missing directories.

calculate_hessian_matrix(atoms, divide_by_masses=False)

Compute the Hessian matrix

The Hessian matrix is the matrix of second derivatives of the potential energy \(\mathcal{V}_\mathrm{int}\) with respect to coordinates, i.e.

\[\frac{\partial^2 \mathcal{V}_\mathrm{int}} {\partial r_{\nu{}i}\partial r_{\mu{}j}},\]

where the indices \(\mu\) and \(\nu\) refer to atoms and the indices \(i\) and \(j\) refer to the components of the position vector \(r_\nu\) along the three spatial directions.

The Hessian matrix has contributions from the pair potential and the embedding energy,

\[\frac{\partial^2 \mathcal{V}_\mathrm{int}}{\partial r_{\nu{}i}\partial r_{\mu{}j}} = \frac{\partial^2 \mathcal{V}_\mathrm{pair}}{ \partial r_{\nu i} \partial r_{\mu j}} + \frac{\partial^2 \mathcal{V}_\mathrm{embed}}{\partial r_{\nu i} \partial r_{\mu j}}.\]

The contribution from the pair potential is

\[\begin{split}\frac{\partial^2 \mathcal{V}_\mathrm{pair}}{ \partial r_{\nu i} \partial r_{\mu j}} &= -\phi_{\nu\mu}'' \left( \frac{r_{\nu\mu i}}{r_{\nu\mu}} \frac{r_{\nu\mu j}}{r_{\nu\mu}} \right) -\frac{\phi_{\nu\mu}'}{r_{\nu\mu}}\left( \delta_{ij}- \frac{r_{\nu\mu i}}{r_{\nu\mu}} \frac{r_{\nu\mu j}}{r_{\nu\mu}} \right) \\ &+\delta_{\nu\mu}\sum_{\gamma\neq\nu}^{N} \phi_{\nu\gamma}'' \left( \frac{r_{\nu\gamma i}}{r_{\nu\gamma}} \frac{r_{\nu\gamma j}}{r_{\nu\gamma}} \right) +\delta_{\nu\mu}\sum_{\gamma\neq\nu}^{N}\frac{\phi_{\nu\gamma}'}{r_{\nu\gamma}}\left( \delta_{ij}- \frac{r_{\nu\gamma i}}{r_{\nu\gamma}} \frac{r_{\nu\gamma j}}{r_{\nu\gamma}} \right).\end{split}\]

The contribution of the embedding energy to the Hessian matrix is a sum of eight terms,

\[\begin{split}\frac{\mathcal{V}_\mathrm{embed}}{\partial r_{\mu j} \partial r_{\nu i}} &= T_1 + T_2 + T_3 + T_4 + T_5 + T_6 + T_7 + T_8 \\ T_1 &= \delta_{\nu\mu}U_\nu'' \sum_{\gamma\neq\nu}^{N}g_{\nu\gamma}'\frac{r_{\nu\gamma i}}{r_{\nu\gamma}} \sum_{\gamma\neq\nu}^{N}g_{\nu\gamma}'\frac{r_{\nu\gamma j}}{r_{\nu\gamma}} \\ T_2 &= -u_\nu''g_{\nu\mu}' \frac{r_{\nu\mu j}}{r_{\nu\mu}} \sum_{\gamma\neq\nu}^{N} G_{\nu\gamma}' \frac{r_{\nu\gamma i}}{r_{\nu\gamma}} \\ T_3 &= +u_\mu''g_{\mu\nu}' \frac{r_{\nu\mu i}}{r_{\nu\mu}} \sum_{\gamma\neq\mu}^{N} G_{\mu\gamma}' \frac{r_{\mu\gamma j}}{r_{\mu\gamma}} \\ T_4 &= -\left(u_\mu'g_{\mu\nu}'' + u_\nu'g_{\nu\mu}''\right) \left( \frac{r_{\nu\mu i}}{r_{\nu\mu}} \frac{r_{\nu\mu j}}{r_{\nu\mu}} \right)\\ T_5 &= \delta_{\nu\mu} \sum_{\gamma\neq\nu}^{N} \left(U_\gamma'g_{\gamma\nu}'' + U_\nu'g_{\nu\gamma}''\right) \left( \frac{r_{\nu\gamma i}}{r_{\nu\gamma}} \frac{r_{\nu\gamma j}}{r_{\nu\gamma}} \right) \\ T_6 &= -\left(U_\mu'g_{\mu\nu}' + U_\nu'g_{\nu\mu}'\right) \frac{1}{r_{\nu\mu}} \left( \delta_{ij}- \frac{r_{\nu\mu i}}{r_{\nu\mu}} \frac{r_{\nu\mu j}}{r_{\nu\mu}} \right) \\ T_7 &= \delta_{\nu\mu} \sum_{\gamma\neq\nu}^{N} \left(U_\gamma'g_{\gamma\nu}' + U_\nu'g_{\nu\gamma}'\right) \frac{1}{r_{\nu\gamma}} \left(\delta_{ij}- \frac{r_{\nu\gamma i}}{r_{\nu\gamma}} \frac{r_{\nu\gamma j}}{r_{\nu\gamma}} \right) \\ T_8 &= \sum_{\substack{\gamma\neq\nu \\ \gamma \neq \mu}}^{N} U_\gamma'' g_{\gamma\nu}'g_{\gamma\mu}' \frac{r_{\gamma\nu i}}{r_{\gamma\nu}} \frac{r_{\gamma\mu j}}{r_{\gamma\mu}}\end{split}\]
Parameters:
  • atoms (ase.Atoms) –

  • divide_by_masses (bool) – Divide block \(\nu\mu\) by \(m_\nu{}m_\mu{}\) to obtain the dynamical matrix

Returns:

D – Block Sparse Row matrix with the nonzero blocks

Return type:

numpy.matrix

Notes

Notation:
  • \(N\) Number of atoms

  • \(\mathcal{V}_\mathrm{int}\) Total potential energy

  • \(\mathcal{V}_\mathrm{pair}\) Pair potential

  • \(\mathcal{V}_\mathrm{embed}\) Embedding energy

  • \(r_{\nu{}i}\) Component \(i\) of the position vector of atom \(\nu\)

  • \(r_{\nu\mu{}i} = r_{\mu{}i}-r_{\nu{}i}\)

  • \(r_{\nu\mu{}}\) Norm of \(r_{\nu\mu{}i}\), i.e.\(\left(r_{\nu\mu{}1}^2+r_{\nu\mu{}2}^2+r_{\nu\mu{}3}^2\right)^{1/2}\)

  • \(\phi_{\nu\mu}(r_{\nu\mu{}})\) Pair potential energy of atoms \(\nu\) and \(\mu\)

  • \(\rho_{\nu}\) Total electron density of atom \(\nu\)

  • \(U_\nu(\rho_nu)\) Embedding energy of atom \(\nu\)

  • \(g_{\delta}\left(r_{\gamma\delta}\right) \equiv g_{\gamma\delta}\) Contribution from atom \(\delta\) to \(\rho_\gamma\)

  • \(m_\nu\) mass of atom \(\nu\)

get_hessian(atoms, format='sparse', divide_by_masses=False)

Calculate the Hessian matrix for a pair potential. For an atomic configuration with N atoms in d dimensions the hessian matrix is a symmetric, hermitian matrix with a shape of (d*N,d*N). The matrix is in general a sparse matrix, which consists of dense blocks of shape (d,d), which are the mixed second derivatives.

Parameters:
  • atoms (ase.Atoms) – Atomic configuration in a local or global minima.

  • format (str, optional) – Output format of the Hessian matrix, either ‘dense’, ‘sparse’ or ‘neighbour-list’. The format ‘sparse’ returns a sparse matrix representations of scipy. The format ‘neighbor-list’ returns a representation within matscipy’s and ASE’s neighbor list format, i.e. the Hessian is returned per neighbor. (Default: ‘dense’)

  • divide_by_masses (bool, optional) – Divided each block entry n the Hessian matrix by sqrt(m_i m_j) where m_i and m_j are the masses of the two atoms for the Hessian matrix.

Returns:

  • If format==’sparse’

  • hessian (scipy.sparse.bsr_matrix) – Hessian matrix in sparse matrix representation

  • If format==’neighbor-list’

  • hessian_ncc (np.ndarray) – Array containing the Hessian blocks per atom pair

  • distances_nc (np.ndarray) – Distance vectors between atom pairs

band_structure()

Create band-structure object for plotting.

calculate_numerical_forces(atoms, d=0.001)

Calculate numerical forces using finite difference.

All atoms will be displaced by +d and -d in all directions.

calculate_numerical_stress(atoms, d=1e-06, voigt=True)

Calculate numerical stress using finite difference.

calculate_properties(atoms, properties)

This method is experimental; currently for internal use.

calculation_required(atoms, properties)
check_state(atoms, tol=1e-15)

Check for any system changes since last calculation.

property directory: str
discard_results_on_any_change = False

Whether we purge the results following any change in the set() method.

export_properties()
get_atoms()
get_birch_coefficients(atoms)

Compute the Birch coefficients (Effective elastic constants at non-zero stress).

Parameters:

atoms (ase.Atoms) – Atomic configuration in a local or global minima.

get_born_elastic_constants(atoms)

Compute the Born elastic constants.

Parameters:

atoms (ase.Atoms) – Atomic configuration in a local or global minima.

get_charges(atoms=None)
get_default_parameters()
get_dipole_moment(atoms=None)
get_dynamical_matrix(atoms)

Compute dynamical matrix (=mass weighted Hessian).

get_elastic_constants(atoms, cg_parameters={'M': None, 'atol': 1e-05, 'callback': None, 'maxiter': None, 'tol': 1e-05, 'x0': None})

Compute the elastic constants at zero temperature. These are sum of the born, the non-affine and the stress contribution.

Parameters:
  • atoms (ase.Atoms) – Atomic configuration in a local or global minima.

  • cg_parameters (dict) –

    Dictonary for the conjugate-gradient solver.

    x0: {array, matrix}

    Starting guess for the solution.

    tol/atol: float, optional

    Tolerances for convergence, norm(residual) <= max(tol*norm(b), atol).

    maxiter: int

    Maximum number of iterations. Iteration will stop after maxiter steps even if the specified tolerance has not been achieved.

    M: {sparse matrix, dense matrix, LinearOperator}

    Preconditioner for A.

    callback: function

    User-supplied function to call after each iteration.

get_forces(atoms=None)
get_magnetic_moment(atoms=None)
get_magnetic_moments(atoms=None)

Calculate magnetic moments projected onto atoms.

get_non_affine_contribution_to_elastic_constants(atoms, eigenvalues=None, eigenvectors=None, pc_parameters=None, cg_parameters={'M': None, 'atol': 1e-05, 'callback': None, 'maxiter': None, 'tol': 1e-05, 'x0': None})

get_non_affine_contribution_to_elastic_constants is deprecated, use elasticity.nonaffine_elastic_contribution instead!

Compute the correction of non-affine displacements to the elasticity tensor. The computation of the occuring inverse of the Hessian matrix is bypassed by using a cg solver.

If eigenvalues and and eigenvectors are given the inverse of the Hessian can be easily computed.

Parameters:
  • atoms (ase.Atoms) – Atomic configuration in a local or global minima.

  • eigenvalues (array) – Eigenvalues in ascending order obtained by diagonalization of Hessian matrix. If given, use eigenvalues and eigenvectors to compute non-affine contribution.

  • eigenvectors (array) – Eigenvectors corresponding to eigenvalues.

  • cg_parameters (dict) –

    Dictonary for the conjugate-gradient solver.

    x0: {array, matrix}

    Starting guess for the solution.

    tol/atol: float, optional

    Tolerances for convergence, norm(residual) <= max(tol*norm(b), atol).

    maxiter: int

    Maximum number of iterations. Iteration will stop after maxiter steps even if the specified tolerance has not been achieved.

    M: {sparse matrix, dense matrix, LinearOperator}

    Preconditioner for A.

    callback: function

    User-supplied function to call after each iteration.

  • pc_parameters (dict) –

    Dictonary for the incomplete LU decomposition of the Hessian.

    A: array_like

    Sparse matrix to factorize.

    drop_tol: float

    Drop tolerance for an incomplete LU decomposition.

    fill_factor: float

    Specifies the fill ratio upper bound.

    drop_rule: str

    Comma-separated string of drop rules to use.

    permc_spec: str

    How to permute the columns of the matrix for sparsity.

    diag_pivot_thresh: float

    Threshold used for a diagonal entry to be an acceptable pivot.

    relax: int

    Expert option for customizing the degree of relaxing supernodes.

    panel_size: int

    Expert option for customizing the panel size.

    options: dict

    Dictionary containing additional expert options to SuperLU.

get_nonaffine_forces(atoms)

Compute the non-affine forces which result from an affine deformation of atoms.

Parameters:

atoms (ase.Atoms) – Atomic configuration in a local or global minima.

get_numerical_non_affine_forces(atoms, d=1e-06)

get_numerical_non_affine_forces is deprecated, use numerical.numerical_nonaffine_forces instead!

Calculate numerical non-affine forces using central finite differences. This is done by deforming the box, rescaling atoms and measure the force.

Parameters:

atoms (ase.Atoms) – Atomic configuration in a local or global minima.

get_potential_energies(atoms=None)
get_potential_energy(atoms=None, force_consistent=False)
get_property(name, atoms=None, allow_calculation=True)

Get the named property.

get_stress(atoms=None)
get_stress_contribution_to_elastic_constants(atoms)

Compute the correction to the elastic constants due to non-zero stress in the configuration. Stress term results from working with the Cauchy stress.

Parameters:

atoms (ase.Atoms) – Atomic configuration in a local or global minima.

get_stresses(atoms=None)

the calculator should return intensive stresses, i.e., such that stresses.sum(axis=0) == stress

ignored_changes: Set[str] = {}

Properties of Atoms which we ignore for the purposes of cache

property label
read(label)

Read atoms, parameters and calculated properties from output file.

Read result from self.label file. Raise ReadError if the file is not there. If the file is corrupted or contains an error message from the calculation, a ReadError should also be raised. In case of succes, these attributes must set:

atoms: Atoms object

The state of the atoms from last calculation.

parameters: Parameters object

The parameter dictionary.

results: dict

Calculated properties like energy and forces.

The FileIOCalculator.read() method will typically read atoms and parameters and get the results dict by calling the read_results() method.

classmethod read_atoms(restart, **kwargs)
reset()

Clear all information from old calculation.

set(**kwargs)

Set parameters like set(key1=value1, key2=value2, …).

A dictionary containing the parameters that have been changed is returned.

Subclasses must implement a set() method that will look at the chaneged parameters and decide if a call to reset() is needed. If the changed parameters are harmless, like a change in verbosity, then there is no need to call reset().

The special keyword ‘parameters’ can be used to read parameters from a file.

set_label(label)

Set label and convert label to directory and prefix.

Examples:

  • label=’abc’: (directory=’.’, prefix=’abc’)

  • label=’dir1/abc’: (directory=’dir1’, prefix=’abc’)

  • label=None: (directory=’.’, prefix=None)

todict(skip_default=True)
property cutoff

matscipy.ios.eam.io module

Read and write tabulated EAM potentials

class matscipy.calculators.eam.io.EAMParameters(symbols, atomic_numbers, atomic_masses, lattice_constants, crystal_structures, number_of_density_grid_points, number_of_distance_grid_points, density_grid_spacing, distance_grid_spacing, cutoff)

Bases: EAMParameters

Embedded Atom Method potential parameters

Parameters:
  • symbols (array_like) – Symbols of the elements coverered by this potential (only for eam/alloy and eam/fs, EMPTY for eam

  • atomic_numbers (array_like) – Atomic numbers of the elements covered by this potential

  • atomic_masses (array_like) – Atomic masses of the elements covered by this potential

  • lattice_constants (array_like) – Lattice constant of a pure crystal with crystal structure as specified in crystal_structures

  • crystal_structures (array_like) – Crystal structure of the pure metal.

  • number_of_density_grid_points (int) – Number of grid points of the embedding energy functional

  • number_of_distance_grid_points (int) – Number of grid points of the electron density function and the pair potential

  • density_grid_spacing (float) – Grid spacing in electron density space

  • distance_grid_spacing (float) – Grid spacing in pair distance space

  • cutoff (float) – Cutoff distance of the potential

atomic_masses

Alias for field number 2

atomic_numbers

Alias for field number 1

count(value, /)

Return number of occurrences of value.

crystal_structures

Alias for field number 4

cutoff

Alias for field number 9

density_grid_spacing

Alias for field number 7

distance_grid_spacing

Alias for field number 8

index(value, start=0, stop=9223372036854775807, /)

Return first index of value.

Raises ValueError if the value is not present.

lattice_constants

Alias for field number 3

number_of_density_grid_points

Alias for field number 5

number_of_distance_grid_points

Alias for field number 6

symbols

Alias for field number 0

matscipy.calculators.eam.io.read_eam(eam_file, kind='eam/alloy')

Read a tabulated EAM potential

There are differnt flavors of EAM, with different storage formats. This function supports a subset of the formats supported by Lammps (http://lammps.sandia.gov/doc/pair_eam.html), * eam (DYNAMO funcfl format) * eam/alloy (DYNAMO setfl format) * eam/fs (DYNAMO setfl format)

Parameters:
  • eam_file (string) – eam alloy file name

  • kind ({'eam', 'eam/alloy', 'eam/fs'}) – kind of EAM file to read

Returns:

  • source (string) – Source informations or comment line for the file header

  • parameters (EAMParameters) – EAM potential parameters

  • F (array_like) – contain the tabulated values of the embedded functions shape = (nb elements, nb of data points)

  • f (array_like) – contain the tabulated values of the density functions shape = (nb elements, nb of data points)

  • rep (array_like) – contain the tabulated values of pair potential shape = (nb elements,nb elements, nb of data points)

matscipy.calculators.eam.io.mix_eam(files, kind, method, f=[], rep_ab=[], alphas=[], betas=[])

mix eam alloy files data set and compute the interspecies pair potential part using the mean geometric value from each pure species

Parameters:
  • files (array of strings) – Contain all the files to merge and mix

  • kind (string) – kinf of eam. Supported eam/alloy, eam/fs

  • method (string, {geometric, arithmetic, weighted, fitted}) – Method used to mix the pair interaction terms. The geometric, arithmetic, and weighted arithmetic average are available. The weighted arithmetic method is using the electron density function values of atom a and b to ponderate the pair potential between species a and \(b\), rep_ab = 0.5(fb/fa * rep_a + fa/fb * rep_b), see [1]. The fitted method is to be used if rep_ab has been previously fitted and is parse as \(rep_ab\) karg.

  • f (np.array) – fitted density term (for FS eam style)

  • rep_ab (np.array) – fitted rep_ab term

  • alphas (array) – fitted alpha values for the fine tuned mixing. rep_ab = alpha_a*rep_a+alpha_b*rep_b

  • betas (array) – fitted values for the fine tuned mixing. f_ab = beta_00*rep_a+beta_01*rep_b f_ba = beta_10*rep_a+beta_11*rep_b

Returns:

  • sources (string) – Source informations or comment line for the file header

  • parameters_mix (EAMParameters) – EAM potential parameters

  • F_ (array_like) – contain the tabulated values of the embedded functions shape = (nb elements, nb elements, nb of data points)

  • f_ (array_like) – contain the tabulated values of the density functions shape = (nb elements, nb elements, nb of data points)

  • rep_ (array_like) – contain the tabulated values of pair potential shape = (nb elements, nb elements, nb of data points)

References

      1. Zhou, R. A. Johnson, and H. N. G. Wadley, Phys. Rev. B, 69, 144113 (2004)

matscipy.calculators.eam.io.write_eam(source, parameters, F, f, rep, out_file, kind='eam')

Write an eam lammps format file

There are differnt flavors of EAM, with different storage formats. This function supports a subset of the formats supported by Lammps (http://lammps.sandia.gov/doc/pair_eam.html), * eam (DYNAMO funcfl format) * eam/alloy (DYNAMO setfl format) * eam/fs (DYNAMO setfl format)

Parameters:
  • source (string) – Source information or comment line for the file header

  • parameters_mix (EAMParameters) – EAM potential parameters

  • F (array_like) – contain the tabulated values of the embedded functions shape = (nb of data points)

  • f (array_like) – contain the tabulated values of the density functions shape = (nb of data points)

  • rep (array_like) – contain the tabulated values of pair potential shape = (nb of data points)

  • out_file (string) – output file name for the eam alloy potential file

  • kind ({'eam', 'eam/alloy', 'eam/fs'}) – kind of EAM file to read

Return type:

None

matscipy.average_atoms.eam.average_atom module

matscipy.calculators.eam.average_atom.average_potential(concentrations, parameters, F, f, rep, kind='eam/alloy', avg_atom='A', atomic_number=999, crystal_structure='unknown', lattice_constant=1.0)

Generate Average-atom potential

The Average-atom (A-atom) potential is a mean-field approximation for random alloys, see Ref. 1. The purpose is to replace the true elements by a single fictious element, the A-atom. A configuration of A-atoms yields approximately the same potential energy as the corresponding random alloy configuration. Other average material properties, e.g. the elastic constants, are reproduced as well. For a full derivation, see Ref. 1.

The A-atom potential has standard EAM form, i.e. it can be tabulated just like any other EAM potential. The potential functions are simply the concentration-weighted averages of the pure element functions:

\[\begin{split}\phi_{AA}\left(r_{\gamma\delta}\right) &= \sum_{X}^{N_T}\sum_{Y}^{N_T} c_{X}c_{Y} \phi_{XY}\left(r_{\gamma\delta}\right) \quad\text{(pair potential)}, \\ U_{A}\left(\rho_\gamma\right) &= \sum_{X}^{N_T}c_{X}U_{X}\left(\rho_\gamma\right) \quad\text{(embedding energy)}, \\ g_A\left(r_{\gamma\delta}\right) &= \sum_{X}^{N_T}c_X g_X\left(r_{\gamma\delta}\right)\quad\text{(electron density)},\;\text{and}\\ m_A &= \sum_{X}^{N_T}c_X m_X\quad\text{(mass)}.\end{split}\]

Note

Currently, only eam/alloy-style potentials can be averaged. The extension to eam/fs should be straightforward, however.

Parameters:
  • concentrations (array_like) – concentrations of the elements in the A-atom

  • parameters (EAMParameters) – EAM potential parameters

  • F (array_like) – tabulated embedding energy functionals

  • f (array_like) – tabulated electron density functions

  • rep (array_like) – tabulated pair potential energy functions

Returns:

  • parameters (EAMParameters) – EAM potential parameters

  • new_F (array_like) – tabulated embedding energy functionals, including A-atom functional

  • new_f (array_like) – tabulated electron density functions, including A-atom function(s)

  • new_rep (array_like) – tabulated pair potential energy functions, including pairs with A-atom

Examples

>>> from matscipy.calculators.eam import io, average_atom
>>> source, parameters, F, f, rep = io.read_eam(
>>>     "ZrCu.onecolumn.eam.alloy"
>>> )
>>> concentrations = [0.5, 0.5]
>>> (new_parameters, new_F, new_f, new_rep) = average_atom.average_potential(
>>>     concentrations, parameters, F, f, rep
>>> )
>>> composition = " ".join(
>>>     [str(c * 100.0) + "% {0},".format(e) for c, e in zip(concentrations, parameters.symbols)]
>>> )
>>> composition = composition.rstrip(",")
>>> source += ", averaged for composition {0}".format(composition)
>>> io.write_eam(
>>>     source,
>>>     new_parameters,
>>>     new_F,
>>>     new_f,
>>>     new_rep,
>>>     "ZrCu.onecolumn.averaged.eam.alloy",
>>>     kind="eam/alloy",
>>> )

Read an EAM potential table for two elements in eam/alloy format, and create a new table with additional A-atom functions for the equicomposition alloy.

References

Notes

Notation:
  • \(N\) Number of atoms

  • \(N_T\) Number of elements

  • \(r_{\nu\mu{}}\) Pair distance of atoms \(\nu\) and \(\mu\)

  • \(\phi_{\nu\mu}(r_{\nu\mu{}})\) Pair potential energy of atoms \(\nu\) and \(\mu\)

  • \(\rho_{\nu}\) Total electron density of atom \(\nu\)

  • \(U_\nu(\rho_nu)\) Embedding energy of atom \(\nu\)

  • \(g_{\delta}\left(r_{\gamma\delta}\right) \equiv g_{\gamma\delta}\) Contribution from atom \(\delta\) to \(\rho_\gamma\)

  • \(m_\nu\) mass of atom \(\nu\)