matscipy.calculators.eam.io
Read and write tabulated EAM potentials
Functions
|
mix eam alloy files data set and compute the interspecies pair potential part using the mean geometric value from each pure species |
|
Read a tabulated EAM potential |
|
Write an eam lammps format file |
Classes
|
Embedded Atom Method potential parameters |
- 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:
EAMParametersEmbedded 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
Methods
count(value, /)Return number of occurrences of value.
index(value[, start, stop])Return first index of value.
- 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=sys.maxsize, /)
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
aandbto ponderate the pair potential between speciesaand \(b\),rep_ab = 0.5(fb/fa * rep_a + fa/fb * rep_b), see [1]. The fitted method is to be used ifrep_abhas 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_bbetas (array) – fitted values for the fine tuned mixing.
f_ab = beta_00*rep_a+beta_01*rep_bf_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
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