Function and class reference¶
Class for handling an atomicrex fitjob. |
|
Class for handling atomic structures. |
This is the Python interface module for atomicrex. It primarily provides convenient access to the atomicrex job object and the atomic structures associated with it.
-
class
atomicrex.
AtomicStructure
¶ Class for handling atomic structures.
-
property
cell
¶ simulation cell metric
- Type
matrix
-
property
cell_origin
¶ origin point of the simulation cell in world coordinates
- Type
vector
-
compute_energy
(compute_forces=False, is_fitting=False, suppress_relaxation=False)¶ Computes the total energy and optionally the forces for this structure.
The function returns a scalar representing the energy of the structure. If the forces are computed they are stored in the structure object.
- Parameters
compute_forces (bool) – Turns on the computation of forces.
is_fitting (bool) – Determines whether this DOF should be relaxed (usually the desired behavior depends on the program phase).
suppress_relaxation (bool) – If True structure relaxation will be suppressed.
-
compute_properties
(is_fitting=False)¶ Computes all enabled properties of the structure.
- Parameters
is_fitting (bool) – Determines whether this DOF should be relaxed (usually the desired behavior depends on the program phase).
-
deactivate_property
(id, fit_enabled=False, output_enabled=False)¶ Deactivates a property in the structure, i.e. removes its contribution to the residual (objective function).
- Parameters
id (string) – property identifier
fit_enabled (bool) – set to True in order to still include the property in the calculation of the objective function
output_enabled (bool) – set to True in order to still include the property when printing the results
-
deform_simulation_cell
(deformation)¶ Applies an affine transformation to the cell and the atoms.
- Parameters
deformation (3x3 matrix) – deformation matrix that is applied to the simulation cell.
-
property
forces
¶ list of atomic forces
- Type
N x 3 doubles
-
get_atoms
(job)¶ Returns the structure as an ASE atoms object.
- Parameters
job – atomicrex job object; this is used to translate the atomicrex atom types to element names as understood by ASE
- Returns
atomic structure
- Return type
ASE atoms object
-
get_cell
()¶ - Returns
cell metric
- Return type
matrix
-
get_forces
()¶ Return the forces on the atoms in the form of an array.
- Returns
atomic forces
- Return type
list of vectors
-
get_number_of_atoms
()¶ - Returns
number of atoms in structure
- Return type
int
-
get_pbc
()¶ Periodic boundary conditions.
- Returns
the periodic boundary conditions
- Return type
list containing three booleans
-
get_positions
()¶ Return the positions of the atoms
- Returns
atomic positions
- Return type
Nx3 array
-
get_potential_energy
()¶ - Returns
potential energy of structure
- Return type
float
-
get_pressure_tensor
(voigtIndex)¶ - Parameters
filename (int) – voigtIndex Index of the pressure tensor component in Voigt notation.
- Returns
pressure tensor entry
- Return type
double
-
get_target_forces
()¶ Get the target forces for a structure.
- Returns
the target forces
- Return type
Nx3 array of floats
-
get_types
()¶ Return the atom types in the form of a list.
Notes
The atom type indices start at 0 (zero).
- Returns
atom types
- Return type
N ints
-
modify_property
(id, target_value=None, relative_weight=1.0, fit_enabled=None, output_enabled=True)¶ Add and/or modify a property.
- Parameters
id (string) – property identifier
target_value (double) – the target value [default: no target value]
relative_weight (double) – relative weight of the property if it is included in the calculation of the objective function
fit_enabled (bool) – calculation of the objective function
output_enabled (bool) – set to True in order to include the property when printing the results
-
property
pbc
¶ Periodic boundary conditions.
- Type
list of three booleans
-
property
positions
¶ atomic positions
- Type
N x 3 doubles
-
property
potential_energy
¶ potential energy of the structure
- Type
float
-
print_properties
(show_all=False)¶ Print properties that are enabled for output and/or fitting.
- Parameters
show_all (bool) – set to True in order to show all properties regardless of their internal preference (as determined from output_enabled and fit_enabled)
Notes
Note that
compute_properties()
must be called first in order to update the stored values of the structure’s properties.
-
property
properties
¶ Properties associated with structure.
- Type
dict
-
set_positions
(positions)¶ Set atomic positions to some new values.
- Parameters
pos (Nx3 array) – new atomic positions
-
set_target_forces
(target_forces)¶ Set the target forces for a structure.
- Parameters
target_forces (Nx3 array of floats) – list of atom vectors representing the target forces
-
set_types
(types)¶ Set atom types.
- Parameters
types (N array of ints) – new atom types
Notes
The atom type indices start at 0 (zero).
atom tags
- Type
list of integers
-
property
target_forces
¶ Get the target forces for a structure.
- Returns
the target forces
- Return type
Nx3 array of floats
-
property
types
¶ atom types
- Type
list of integers
-
write_lammps_dump
(filename, ghosts=False)¶ Exports the structure to a LAMMPS dump file.
- Parameters
filename (string) – Name of the output file.
ghosts (bool) – If True, ghost atoms used during computation of the energy will be included in the output.
-
write_poscar
(filename, ghosts=False)¶ Exports the structure to a POSCAR file.
- Parameters
filename (string) – Name of the output file.
ghosts (bool) – If True, ghost atoms used during computation of the energy will be included in the output.
-
property
-
class
atomicrex.
Job
¶ Class for handling an atomicrex fitjob.
-
add_ase_structure
(structure_id, atoms, structure_group=None)¶ Add a structure based on an ASE atoms object.
- Parameters
structure_id (string) – structure name; to be used later for accessing structure in dictionary
atoms (ASE atoms object) – atomic configuration
structure_group (reference to a structure group) – reference to structure group, to which this structure is to be assigned; the default (None) implies the root group
-
add_library_structure
(id, structure, params, structure_group=None)¶ Add a structure from the library of predefined structures.
- Parameters
id (string) – structure name that will be used later to access the structure in the dictionary
structure (string) –
type of structure
single element: FCC, BCC, DIA, SC, HCP, DHCP, beta tin, omega etc.
multiple elements: rocksalt, cesium chloride, zinc blende, D8a, L12, fluorite, C15, Bh, D2d, Ni17Th2, Th2Zn17, L10, wurtzite etc.
params (dictionary) – structural parameters including e.g, alat, clat, ca_ratio, type, type_A, type_B etc.
structure_group (reference to a structure group) – reference to structure group, to which this structure is to be assigned; the default (None) implies the root group
-
calculate_gradient
(parameters=None, eps=0.0001)¶ Computes the gradient of the objective function with respect to the potential parameters using a central finite difference scheme.
- Parameters
parameters (list of floats) – the parameter set, at which the derivative is to be computed
eps – step length for numerical differentiation
-
calculate_hessian
(parameters=None, eps=0.0001)¶ Computes the matrix of second derivative of the objective function with respect to the potential parameters using a central finite difference scheme.
- Parameters
parameters (list of floats) – the parameter set, at which the derivative is to be computed
eps – step length for numerical differentiation
-
calculate_residual
(params=None, style='squared')¶ Calculates the the objective function to be minimized during fitting. This is the weighted sum of the deviations of all properties from their target values.
- Parameters
params (list) – parameter vector, if specified it will overwrite the current parameter vector associated with this job
style (string) – method employed for calculating the residual contribution of this property to the objective function; possible values are: squared, squared_relative, user_defined, absolute_diff
-
property
derived_properties
¶ derived-properties associated with job.
- Type
dict
-
get_potential_lower_bounds
()¶ - Returns
potential parameter lower bounds (DOFs) associated with job
- Return type
list of floats
-
get_potential_parameters
()¶ - Returns
potential parameters (DOFs) associated with job
- Return type
list of floats
-
get_potential_upper_bounds
()¶ - Returns
potential parameter upper bounds (DOFs) associated with job
- Return type
list of floats
-
property
name
¶ name of fitjob
- Type
str
-
output_results
()¶ Write resulting potential parameters to file.
-
parse_input_file
(filename)¶ Parse the general and potential parameters from an XML file.
- Parameters
filename (str) – name of XML input file
-
perform_fitting
()¶ Perform the actual fitting by minimizing the residual.
-
property
potentials
¶ potentials associated with job.
- Type
dict
-
prepare_fitting
()¶ Set up a list of all degrees of freedom and properties that are to be fitted. This function needs to be called prior to any computation. It is called automatically before a model is fitted via the
perform_fitting()
call.
-
print_potential_parameters
(only_active=False)¶ Print potential parameters associated with job. By default all parameters (including inactive ones) are written to stdout.
- Parameters
only_active (bool) – If True exclude inactive parameters
-
print_properties
()¶ Print properties associated with job.
-
set_potential_parameters
(values)¶ - Parameters
the potential parameters (DOFs) associated with the job. (Sets) –
values (array) – list of double values representing the parameters
-
set_verbosity
(verbosity)¶ Set the verbosity.
- Parameters
verbosity (int) – verbosity level (0=none, 1=minimum, 2=medium, 3=maximum, 4=debug)
-
property
structures
¶ structures associated with job.
- Type
dict
-