matscipy.cauchy_born

Classes

CubicCauchyBorn(el, a0, calc[, lattice])

Corrector for the cauchy-born prediction of atomistic positions in multi-lattices subject to continuum strain This model exploits the symmetry triad down the 111 axis that all cubic crystals possess by making models that only fit a single component of the shift-correction vector.

RegressionModel()

Simple regression model object wrapper for np.linalg.lstsq for predicting shifts.

class matscipy.cauchy_born.RegressionModel

Bases: object

Simple regression model object wrapper for np.linalg.lstsq for predicting shifts.

Methods

fit(phi, values)

Fit a simple least-squares regression model to data.

load()

Loads the regression model from file: CB_model.txt

predict(phi)

Evaluate the simple least-squares regression model.

predict_gradient(grad_phi)

Evaluate the simple least-squares regression model using the derivative of the basis functions to get the gradient.

save()

Saves the regression model to file: CB_model.txt

fit(phi, values)

Fit a simple least-squares regression model to data. Saves output vector to self.model.

Parameters:
  • phi (array_like, 2D) – Design matrix of shape [number of data points, number of basis functions]. Contains each basis function evaluated at each data point.

  • values (array_like) – value of function to be fitted evaluated at each data point

predict(phi)

Evaluate the simple least-squares regression model.

Parameters:

phi (array_like, 2D) – Design matrix of shape [number of data points, number of basis functions]. Contains each basis function evaluated at each data point.

Returns:

predictions – 1D vector of least squares model predictions for simple model.

Return type:

array_like

predict_gradient(grad_phi)

Evaluate the simple least-squares regression model using the derivative of the basis functions to get the gradient.

Parameters:

grad_phi (array_like, 2D) – Design matrix of shape [number of data points, number of basis functions]. Contains each gradient basis function evaluated at each data point.

Returns:

predictions – 1D vector of least squares model gradient predictions for simple model.

Return type:

array_like

save()

Saves the regression model to file: CB_model.txt

load()

Loads the regression model from file: CB_model.txt

class matscipy.cauchy_born.CubicCauchyBorn(el, a0, calc, lattice=<ase.lattice.cubic.DiamondFactory object>)

Bases: object

Corrector for the cauchy-born prediction of atomistic positions in multi-lattices subject to continuum strain This model exploits the symmetry triad down the 111 axis that all cubic crystals possess by making models that only fit a single component of the shift-correction vector.

The other 2 components can be obtained by calling the same model with strain states that are rotated about the 111 axis.

Methods

apply_shifts(atoms, shifts[, mask])

Apply predicted Cauchy-Born corrector shifts to atoms object with an optional mask.

basis_function_evaluation(E_vecs)

Function that evaluates the polynomial basis functions used in the regression model on a set of strain vectors to build the design matrix phi.

check_for_refit(A, atoms, forces[, tol, ...])

Check and perform a refit of the Cauchy-Born shift regression model based on the error metric defined above, adding additional strain-states to the dataset from the highest force atoms.

eval_shift(E_vec, unitcell)

Function that calculates the cauchy-born shift corrector for a given Green-Lagrange strain tensor.

evaluate_E(A, atoms, E_func[, cell, coordinates])

Get the Green-Lagrange strain tensor field for a system of atoms in the lab frame directly from a provided function.

evaluate_F(A, atoms, F_func[, cell, coordinates])

Get the deformation gradient tensor field for a system of atoms in the lab frame.

evaluate_F_or_E(A, atoms[, F_func, E_func, ...])

Get the deformation gradient tensor field or Green-Lagrange strain tensor field for a system of atoms in the lab frame, depending on which function is provided.

evaluate_shift_gradient_regression(E, dE)

Evaluate the gradient of the regression model a given set of Green-Lagrange strain tensors and their gradients.

evaluate_shift_model(E[, method])

Evaluate a fitted shift model and return the shift predictions for each atom.

fit_taylor([de])

Fit a simple taylor-expansion type model to predict the cauchy-born shift corrector at a given applied strain using finite differences.

get_cb_error(atoms[, mask, forces])

Get an approximate quantifier of the error in an atoms object due to poor or not implemented Cauchy-Born corrector shifts.

get_data_points(E_vec)

Get the 'true' shifts (by building and relaxing a unit cell under that strain state) using a series of supplied strain vectors.

get_shift_gradients(A, atoms[, de, F_func, ...])

Get the deformation gradient tensor field or Green-Lagrange strain tensor field for a system of atoms in the lab frame, depending on which function is provided.

grad_basis_function_evaluation(E_vecs, dE_vecs)

Function that evaluates the derivative of the polynomial basis functions used in the regression model on a set of strain vectors to build the gradient design matrix grad_phi.

initial_regression_fit([initial_samples])

Perform an initial fit of the regression model with some number of initial samples.

load_regression_model()

Loads the regression model from file: CB_model.txt

load_taylor()

Read the results of the taylor expansion rather than re-calculating each time.

permutation(strain_vec, perm_shift)

Rotate a Green-Lagrange strain tensor in Voigt notation about the 111 axis a number of times given by perm_shift.

predict_shifts(A, atoms[, F_func, E_func, ...])

Get the deformation gradient tensor field or Green-Lagrange strain tensor field for a system of atoms in the lab frame, depending on which function is provided.

refit_regression(atoms, E_voigt)

Refit the regression model to a new set of data.

save_regression_model()

Saves the regression model to file: CB_model.txt

save_taylor()

Save the results of the taylor expansion in to numpy readable files rather than re-calculating each time.

set_sublattices(atoms, A[, read_from_atoms])

Apply a small strain to all atoms in the supplied atoms structure and determine which atoms belong to which sublattice using forces.

switch_sublattices(atoms)

Switch the sublattice masks on a set of atoms around such that lattice1mask = lattice2mask and vice versa

tensor_field_3D_from_atoms(atoms, func[, ...])

Generate a 3D tensor field (F or E) for an atoms object in a given coordinate system, using a provided function

get_model

set_model

__init__(el, a0, calc, lattice=<ase.lattice.cubic.DiamondFactory object>)
Parameters:
  • el (string) – ASE chemical element that constitutes lattice

  • a0 (float) – Lattice constant for cubic lattice

  • calc (ASE calculator object) – ASE calculator.

  • lattice (ASE lattice builder function) – ASE function from ase.lattice.cubic

set_sublattices(atoms, A, read_from_atoms=False)

Apply a small strain to all atoms in the supplied atoms structure and determine which atoms belong to which sublattice using forces. NOTE as this method is based on forces, it does not work in cases where atom forces are high already, This means it should be applied to an unstrained version of the atoms structure and any sublattices assigned adjacent to free surfaces cannot be trusted. This function updates self.lattice1mask and self.lattice2mask with the measured sublattices.

Parameters:
  • atoms (ASE atoms object) – atoms structure with atoms on seperate sublattices that need assigning

  • A (3x3 numpy array, floats) – rotation matrix of the form [x^T,y^T,z^T], where x^T,y^T,z^T are normalised column vectors of the lab frame directions expressed in terms of the crystal lattice directions.

  • read_from_atoms (bool) – whether or not to read the sublattices directly from the atoms object provided, rather than working them out by applying a small strain and measuring forces. The user will have to set the sublattices themselves when building the atoms structure

switch_sublattices(atoms)

Switch the sublattice masks on a set of atoms around such that lattice1mask = lattice2mask and vice versa

Parameters:

atoms (ASE atoms object with defined sublattices) –

fit_taylor(de=0.0001)

Fit a simple taylor-expansion type model to predict the cauchy-born shift corrector at a given applied strain using finite differences. Sets self.grad_f to the vector of first derivatives of the first shift component with each of the 6 strain components and self.hess_f to the hessian of the first shift component with the strain components. NOTE this model is generally worse than the model provided later in initial_regression_fit.

Parameters:

de (float) – tolerance for finite difference gradient calculations

save_taylor()

Save the results of the taylor expansion in to numpy readable files rather than re-calculating each time.

load_taylor()

Read the results of the taylor expansion rather than re-calculating each time.

tensor_field_3D_from_atoms(atoms, func, coordinates='cart3D', cell=None, *args, **kwargs)

Generate a 3D tensor field (F or E) for an atoms object in a given coordinate system, using a provided function

Parameters:
  • atoms (ASE atoms object) – Atoms object where atoms sit at the coordinates needed to generate the tensor field from func.

  • func (Function) – Function to calculate the tensor field from a set of coordinates. The system of coordinates that the function accepts is given by ‘coordinates’ and the function also accepts anything extra specified in *args and **kwargs. Function must return the field for the atoms in form of a numpy array shape [natoms,ndims,ndims], where output is specified in cartesian coordinates.

  • coordinates (string) – The coordinates system that func accepts. Must be ‘cylind2D’, ‘cylind3D’, ‘spherical’, ‘cart2D’, ‘cart3D’

  • cell (array_like, optional) – An optional argument allowing the user to specify the cell parameters used to define the coordinate system.

Returns:

field3D – The 3D tensor field for all the atoms in the system, expressed in the lab frame.

Return type:

array_like

evaluate_F_or_E(A, atoms, F_func=None, E_func=None, cell=None, coordinates='cart3D', *args, **kwargs)

Get the deformation gradient tensor field or Green-Lagrange strain tensor field for a system of atoms in the lab frame, depending on which function is provided. From these, find the Green-Lagrange strain tensor field in the lattice frame.

Parameters:
  • A (3x3 numpy array, floats) – rotation matrix of the form [x^T,y^T,z^T], where x^T,y^T,z^T are normalised column vectors of the lab frame directions expressed in terms of the crystal lattice directions.

  • atoms (ASE atoms object) – Atoms object where atoms sit at the coordinates needed to generate the deformation gradient tensor field from F_func.

  • F_func (Function) – Function to calculate the deformation gradient tensor field from a set of atomic coordinates. The system of coordinates that the function accepts is given by ‘coordinates’ and the function also accepts anything extra specified in *args and **kwargs. Function must return the deformation gradient field for the atoms in form of a numpy array shape [natoms,ndims,ndims], where F is specified in cartesian coordinates. The returned tensor field must be in the lab frame Optional, but one of F_func or E_func must be provided.

  • E_func (Function) – Function to calculate the Green-Lagrange tensor field from a set of atomic coordinates. The system of coordinates that the function accepts is given by ‘coordinates’ and the function also accepts anything extra specified in *args and **kwargs. Function must return the Green-Lagrange tensor field for the atoms in form of a numpy array shape [natoms,ndims,ndims], where E is specified in cartesian coordinates, and the returned tensor field is in the lab frame. Optional, but one of F_func or E_func must be provided.

  • cell (array_like, optional) – An optional argument allowing the user to specify the cell dimensions which will be used to find the centre of the cell. If not provided, the centre of the cell will be found from the atoms object provided.

  • coordinates (string) – Specifies the coordinate system that the function (F_func or E_func) provided accepts. The coordinates of the atoms in atoms will be converted to this coordinate system before being passed to the function as a vector. Default is 3D cartesian coordinates. Other options are cart2D, cart3D, cylind2D, cylind3D, spherical

Returns:

evaluate_F(A, atoms, F_func, cell=None, coordinates='cart3D', *args, **kwargs)

Get the deformation gradient tensor field for a system of atoms in the lab frame. From this, find the Green-Lagrange strain tensor field in the lattice frame.

Parameters:
  • A (3x3 numpy array, floats) – rotation matrix of the form [x^T,y^T,z^T], where x^T,y^T,z^T are normalised column vectors of the lab frame directions expressed in terms of the crystal lattice directions.

  • atoms (ASE atoms object) – Atoms object where atoms sit at the coordinates needed to generate the deformation gradient tensor field from F_func.

  • F_func (Function) – Function to calculate the deformation gradient tensor field from a set of atomic coordinates. The system of coordinates that the function accepts is given by ‘coordinates’ and the function also accepts anything extra specified in *args and **kwargs. Function must return the deformation gradient field for the atoms in form of a numpy array shape [natoms,ndims,ndims], where F is specified in cartesian coordinates. The returned tensor field must be in the lab frame.

  • cell (array_like, optional) – An optional argument allowing the user to specify the cell dimensions which will be used to find the centre of the cell. If not provided, the centre of the cell will be found from the atoms object provided.

  • coordinates (string) – Specifies the coordinate system that the function (F_func or E_func) provided accepts. The coordinates of the atoms in atoms will be converted to this coordinate system before being passed to the function as a vector. Default is 3D cartesian coordinates. Other options are cart2D, cart3D, cylind2D, cylind3D, spherical

Returns:

evaluate_E(A, atoms, E_func, cell=None, coordinates='cart3D', *args, **kwargs)

Get the Green-Lagrange strain tensor field for a system of atoms in the lab frame directly from a provided function.

Parameters:
  • A (3x3 numpy array, floats) – rotation matrix of the form [x^T,y^T,z^T], where x^T,y^T,z^T are normalised column vectors of the lab frame directions expressed in terms of the crystal lattice directions.

  • atoms (ASE atoms object) – Atoms object where atoms sit at the coordinates needed to generate the deformation gradient tensor field from F_func.

  • E_func (Function) – Function to calculate the Green-Lagrange tensor field from a set of atomic coordinates. The system of coordinates that the function accepts is given by ‘coordinates’ and the function also accepts anything extra specified in *args and **kwargs. Function must return the Green-Lagrange tensor field for the atoms in form of a numpy array shape [natoms,ndims,ndims], where E is specified in cartesian coordinates, and the returned tensor field is in the lab frame.

  • cell (array_like, optional) – An optional argument allowing the user to specify the cell dimensions which will be used to find the centre of the cell. If not provided, the centre of the cell will be found from the atoms object provided.

  • coordinates (string) – Specifies the coordinate system that the function (F_func or E_func) provided accepts. The coordinates of the atoms in atoms will be converted to this coordinate system before being passed to the function as a vector. Default is 3D cartesian coordinates. Other options are cart2D, cart3D, cylind2D, cylind3D, spherical

Returns:

predict_shifts(A, atoms, F_func=None, E_func=None, coordinates='cart3D', method='regression', *args, **kwargs)

Get the deformation gradient tensor field or Green-Lagrange strain tensor field for a system of atoms in the lab frame, depending on which function is provided. From these, find the Green-Lagrange strain tensor field in the lattice frame.

Parameters:
  • A (3x3 numpy array, floats) – rotation matrix of the form [x^T,y^T,z^T], where x^T,y^T,z^T are normalised column vectors of the lab frame directions expressed in terms of the crystal lattice directions.

  • atoms (ASE atoms object) – Atoms object where atoms sit at the coordinates needed to generate the deformation gradient tensor field from F_func.

  • F_func (Function) – Function to calculate the deformation gradient tensor field from a set of atomic coordinates. The system of coordinates that the function accepts is given by ‘coordinates’ and the function also accepts anything extra specified in *args and **kwargs. Function must return the deformation gradient field for the atoms in form of a numpy array shape [natoms,ndims,ndims], where F is specified in cartesian coordinates. The returned tensor field must be in the lab frame Optional, but one of F_func or E_func must be provided.

  • E_func (Function) –

    Function to calculate the Green-Lagrange tensor field from a set of atomic coordinates. The system of coordinates that the function accepts is given by ‘coordinates’ and the function also accepts anything extra specified in *args and **kwargs. Function must return the Green-Lagrange tensor field for the atoms in form of a numpy array shape [natoms,ndims,ndims], where E is specified in cartesian coordinates, and the returned tensor field is in the lab frame.

    Optional, but one of F_func or E_func must be provided.

  • coordinates (string) – Specifies the coordinate system that the function (F_func or E_func) provided accepts. The coordinates of the atoms in atoms will be converted to this coordinate system before being passed to the function as a vector. Default is 3D cartesian coordinates. Other options are cart2D, cart3D, cylind2D, cylind3D, spherical

  • method (string) – Specifies which method to use when evaluating shifts. Options are ‘taylor’ (for the taylor expansion model) and ‘regression’ for the regression model. ‘regression’ by default.

Returns:

shifts – The cauchy-born corrector shift predictions made by the model

Return type:

array-like

permutation(strain_vec, perm_shift)

Rotate a Green-Lagrange strain tensor in Voigt notation about the 111 axis a number of times given by perm_shift. This corresponds to permuting the top three components and bottom three components seperately. For example, the vector [1,2,3,4,5,6] becomes [2,3,1,5,6,4] under this operation. Due to the triad 111 symmetry in diamond, feeding these rotated strain states to the shift prediction model allows the prediction of the full shift vector.

Parameters:
  • strain_vec (array_like) – 2D array of size [natoms,6], where each row is a Green-Lagrange strain tensor written in Voigt notation.

  • perm_shift (integer) – The permeutation number wanted, with 1 corresponding to a 120 degree rotation about 111 (these strain vectors give the y component of shift prediction when fed to the model) and 2 corresponding to a 240 degree rotation about 111 (these strain vectors give the z component of shift prediction when fed to the model)

Returns:

strain_vec_perm – Rotated version of strain_vec

Return type:

array_like

evaluate_shift_model(E, method='regression')

Evaluate a fitted shift model and return the shift predictions for each atom.

Parameters:
  • E (array_like) – list of Green-Lagrange strain tensors defining the strain state on each atom, of the form [natoms,3,3]

  • method (string) – Specifies which method to use when evaluating shifts. Options are ‘taylor’ (for the taylor expansion model) and ‘regression’ for the regression model. ‘regression’ by default.

Returns:

shifts – The cauchy-born corrector shift predictions made by the model specified in ‘method’

Return type:

array-like

apply_shifts(atoms, shifts, mask=None)

Apply predicted Cauchy-Born corrector shifts to atoms object with an optional mask.

Parameters:
  • atoms (ASE atoms object) – Atoms object to apply Cauchy-Born corrector shifts too.

  • shifts (array_like) – Cauchy-Born corrector shifts to apply

  • mask (array_like, bool, optional) – Optional mask which specifies which atoms to apply shifts to. All atoms are shifted by default.

get_cb_error(atoms, mask=None, forces=None)

Get an approximate quantifier of the error in an atoms object due to poor or not implemented Cauchy-Born corrector shifts.

Parameters:
  • atoms (ASE atoms object) – ASE atoms for which to quantify the Cauchy-Born corrector related error.

  • mask (array_like, bool) – Mask for which atoms to include in the error quantification process.

  • forces (array_like) – Forces on atoms. Can be provided in order to stop this function from carrying out a potentially expensive atoms.get_forces() operation.

Returns:

predictions – A quantification of the Cauchy-Born shift related error in an atoms object, which can be used as an objective to minimise.

Return type:

cb_err

initial_regression_fit(initial_samples=10)

Perform an initial fit of the regression model with some number of initial samples.

Parameters:

initial_samples (int) – The number of Latin Hypercube samples to fit the initial Cauchy-Born regression model to.

check_for_refit(A, atoms, forces, tol=0.001, mask=None, F_func=None, E_func=None, coordinates='cart3D', refit_points=10, err_vec=None, *args, **kwargs)

Check and perform a refit of the Cauchy-Born shift regression model based on the error metric defined above, adding additional strain-states to the dataset from the highest force atoms.

Parameters:
  • A (3x3 numpy array, floats) – Rotation matrix of the form [x^T,y^T,z^T], where x^T, y^T, z^T are normalised column vectors of the lab frame directions expressed in terms of the crystal lattice directions.

  • atoms (ASE atoms object) – Atoms object where atoms sit at the coordinates needed to generate the deformation gradient tensor field from F_func.

  • forces (array_like) – The forces on the atoms.

  • tol (float) – The acceptable level of Cauchy-Born error per atom.

  • mask (array_like, bool) – Mask for which atoms to include in the error quantification process.

  • F_func (Function) – Function to calculate the deformation gradient tensor field from a set of atomic coordinates. The system of coordinates that the function accepts is given by ‘coordinates’ and the function also accepts anything extra specified in *args and **kwargs. Function must return the deformation gradient field for the atoms in form of a numpy array shape [natoms,ndims,ndims], where F is specified in cartesian coordinates. The returned tensor field must be in the lab frame. Optional, but one of F_func or E_func must be provided.

  • E_func (Function) – Function to calculate the Green-Lagrange tensor field from a set of atomic coordinates. The system of coordinates that the function accepts is given by ‘coordinates’ and the function also accepts anything extra specified in *args and **kwargs. Function must return the Green-Lagrange tensor field for the atoms in form of a numpy array shape [natoms,ndims,ndims], where E is specified in cartesian coordinates, and the returned tensor field is in the lab frame. Optional, but one of F_func or E_func must be provided.

  • coordinates (string) – Specifies the coordinate system that the function (F_func or E_func) provided accepts. The coordinates of the atoms in atoms will be converted to this coordinate system before being passed to the function as a vector. Default is 3D cartesian coordinates. Other options are cart2D, cart3D, cylind2D, cylind3D, spherical

  • refit_points (int) – The number of strain states from high force atoms to add to the dataset for refitting.

  • err_vec (array_like) – A vector of the errors on each atom which can be used in place of the absolute force on each atom as a criteria for selecting the 10 highest error strain states to re-fit the Cauchy-Born model to.

Returns:

  • refit (integer, 1 or 0) – 0 if no refit took place, 1 if it did.

  • cb_err_per_atom (float) – The computed Cauchy-Born shift error metric per atom

refit_regression(atoms, E_voigt)

Refit the regression model to a new set of data.

Parameters:
  • atoms (ASE atoms object) – Atoms object where atoms sit at the coordinates needed to generate the deformation gradient tensor field from F_func.

  • E_voigt (array_like) – Numpy array of the Green-Lagrange strain tensors (in Voigt notation) to refit to.

get_data_points(E_vec)

Get the ‘true’ shifts (by building and relaxing a unit cell under that strain state) using a series of supplied strain vectors. The cauchy-born shift correction along all 3 axes can be measured from one applied strain state and relaxation. The regression model only predicts the component of the shift vector parallel to axis 1, but we can still use all 3 predicted values for fitting if we find the rotated strain tensors about 111 which would predict the shifts along axis 2 and 3.

Parameters:

E_vec (array_like) – Numpy array of the Green-Lagrange strain tensors (in Voigt notation) to apply

Returns:

  • strains (array_like) – A larger array including the supplied E_vec array, as well as the versions of E_vec which have been rotated about 120 degrees about the [111] axis.

  • shift_vals (array_like) – An array of shift components along the 1 axis corresponding to the matching strain state in ‘strains’.

eval_shift(E_vec, unitcell)

Function that calculates the cauchy-born shift corrector for a given Green-Lagrange strain tensor.

Parameters:
  • E_vec (array_like) – Green-Lagrange strain state to apply in voigt notation

  • unitcell (ASE atoms object) – Single unit cell of crystal to apply strain vector to and to find shift

Returns:

shift_diff – Vector of cauchy-born shift correction for applied strain state

Return type:

array_like

basis_function_evaluation(E_vecs)

Function that evaluates the polynomial basis functions used in the regression model on a set of strain vectors to build the design matrix phi.

Parameters:

E_vecs (array_like) – Set of Green-Lagrange strain states in voigt notation to fit to

Returns:

phi – Design matrix composed of all of the basis functions evaluated at every data point given in E_vecs.

Return type:

array_like

grad_basis_function_evaluation(E_vecs, dE_vecs)

Function that evaluates the derivative of the polynomial basis functions used in the regression model on a set of strain vectors to build the gradient design matrix grad_phi.

Parameters:
  • E_vecs (array_like) – Set of Green-Lagrange strain states in voigt notation to fit to

  • dE_vecs (array_like) – Set of vectors holding the derivatives of E_vecs with respect to some other parameter

Returns:

grad_phi – Design matrix composed of all of the gradient basis functions evaluated at every data point given in E_vecs.

Return type:

array_like

evaluate_shift_gradient_regression(E, dE)

Evaluate the gradient of the regression model a given set of Green-Lagrange strain tensors and their gradients.

Parameters:
  • E (array_like) – 2D array of size [natoms,6], where each row is a Green-Lagrange strain tensor written in Voigt notation.

  • dE (array_like) – 2D array of size [natoms, 6], where each row is the derivative of the corresponding row in E with respect to some parameter

Returns:

predictions – prediction of gradient of the Cauchy-Born shift corrector model at each of the given E points.

Return type:

array_like

get_shift_gradients(A, atoms, de=1e-05, F_func=None, E_func=None, coordinates='cart3D', *args, **kwargs)

Get the deformation gradient tensor field or Green-Lagrange strain tensor field for a system of atoms in the lab frame, depending on which function is provided. From these, get the gradients of the Cauchy-Born shift corrector field field in the lattice frame. The supplied F_func or E_func needs to take parameter de, used for finite differences to evaluate the derivative of the given field

Parameters:
  • A (3x3 numpy array, floats) – rotation matrix of the form [x^T,y^T,z^T], where x^T,y^T,z^T are normalised column vectors of the lab frame directions expressed in terms of the crystal lattice directions.

  • atoms (ASE atoms object) – Atoms object where atoms sit at the coordinates needed to generate the deformation gradient tensor field from F_func.

  • F_func (Function) – Function to calculate the deformation gradient tensor field from a set of atomic coordinates. The system of coordinates that the function accepts is given by ‘coordinates’ and the function also accepts anything extra specified in *args and **kwargs. Function must return the deformation gradient field for the atoms in form of a numpy array shape [natoms,ndims,ndims], where F is specified in cartesian coordinates. The returned tensor field must be in the lab frame. Optional, but one of F_func or E_func must be provided. The function must accept de as a keyword argument for calculating the gradient of the field by finite differences.

  • E_func (Function) – Function to calculate the Green-Lagrange tensor field from a set of atomic coordinates. The system of coordinates that the function accepts is given by ‘coordinates’ and the function also accepts anything extra specified in *args and **kwargs. Function must return the Green-Lagrange tensor field for the atoms in form of a numpy array shape [natoms,ndims,ndims], where E is specified in cartesian coordinates, and the returned tensor field is in the lab frame. Optional, but one of F_func or E_func must be provided. The function must accept de as a keyword argument for calculating the gradient of the field by finite differences.

  • coordinates (string) – Specifies the coordinate system that the function (F_func or E_func) provided accepts. The coordinates of the atoms in atoms will be converted to this coordinate system before being passed to the function as a vector. Default is 3D cartesian coordinates. Other options are cart2D, cart3D, cylind2D, cylind3D, spherical

Returns:

shifts – The cauchy-born corrector shift predictions made by the model

Return type:

array-like

save_regression_model()

Saves the regression model to file: CB_model.txt

load_regression_model()

Loads the regression model from file: CB_model.txt

get_model()
set_model(model)