principia_materia.phonon_id.frozen_phonons module

class principia_materia.phonon_id.frozen_phonons.FrozenPhonons(structure, Qpoint, supa=None, hidden_order=1, fdtype='c', tol=1e-06)

Bases: object

Compute coefficients of phonons and their interactions.

The abstract class with the core functionalities to compute phonons and their interactions. For now the symmetry is not considered, it will be in the LID variation, since lone-Q approach can be used to compute the interaction coefficients without symmetry.

Parameters:
  • structure (CrystalFTG object or a structure file of CrystalFTG) – The crystal structure.

  • Qpoint (array of Fraction or a parse_array supported Fraction array.) – The Q-point

  • supa (array of int, shape(dim, dim)) – The supercell matrix.

  • hidden_order (int, optional, default to 1) – The order of derivatives from first principles.

  • fdtype (str, choice of ["c", "f", "b"]) – Type of finite displacement.

  • tol (float, optional, default to 1.0E-6) – Error tolerence.

check_supercell()

Check if supercell can accommodates the Q-point.

create_jobs(job_handler, deltas, prefix=None, delta_format=<function format_delta_dirname>, measurement_label=<function format_measurement_label>, jobname_prefix='', config_file='finite_displacements.yml', fdseries_filename='fdseries.hdf5', append=False, dry_run=False)

Create compute jobs for the frozen phonon calculations.

Each measurement has its own FDSeries class to allow flexibility for different measurements to have different deltas, and each measurement can be managed independently to create extra calculations for better quality.

Parameters:
  • job_handler (ComputeJobSeries or JobsDB) – Either a ComputeJobSeries or a JobsDB object as the interface to the jobs.

  • deltas (list of numbers, list of arrays, 1-d or 2-d array) –

    Displacement amplitude. Can take following types:

    1. a number, will be broadcast into a shape(n_measurements, 1) array.

    2. an 1-d array of length N, will be broadcast into a shape(n_measurements, N) array.

    3. a 2-d array of shape(n_measurements, N).

  • prefix (str, optional, default to None) – Prefix to the jobnames.

  • delta_format (callable, optional, default to function: format_delta_dirname.) – A function that format a delta into a string for directory name. Must be in form of def function(delta, ...): ....

  • jobname_prefix (str) – Prefix for jobnames.

  • config_file (str, optional, default to "finite_displacements.yml") – The name of the finite displacements configuration file.

  • fdseries_filename (str, optional, default to "fdseries.hdf5") – The name of the FDSeries data file to store finite displacement information.

  • dry_run (bool, optional, default to False) – If True, don’t actually create the job.

property derivative_order

The order of the finite displacements derivative.

get_dynamic_tensor()
property order

The order of the phonon interaction.

set_displacements()

Abstract method to construct displacements vectors for finite displacements.

set_displacements_basis()

Abstract method to construct the displacements basis vectors.

set_errortail_results(pick_min=3, pick_max=None, consecutive=False, penalty=<function penalty_linear_mse>, separate_complex=True, output=None, overwrite=True)

Compute error tail for the FrozenPhonons result.

Parameters:
  • pick_min (int, optional, default to 3) – The minimum number of picks for delta selection.

  • pick_max (int, optional, default to None) – The maximum number of picks for delta selection. If None, the pick n scheme is used with n=pick_min. If not None, the pick N schcme is used with N=[pick_min, pick_max].

  • consecutive (bool, optional, default to False) – Whether to pick consecutive deltas in the picking process.

  • penalty (callable, optional, default to penalty_linear_mse) – The penalty function to determine the best fit.

  • separate_complex (bool, optional, default to False) – Whether to fit real and imaginary part of the complex data points separetely.

  • output (str, optional, default to None) – Path of the output file. If not None, save the output of errortail fit to the file.

_fp_errortail

The errortail extrapolated finite displacements results of the FrozenPhonons calculations.

Type:

array of complex

set_realspace_displacements(allreal=False)

Find real space displacements from the complex reciprocal space displacements.

realdisp_cs_index

The list of cos/sin pairs (real/imag pairs) that are required to compute the complex displacements of the measurement.

Type:

list of arrays of int

realdisp_transformation_matrices

The transformation matrices from the complex displacement derivatives to the real space displacement derivatives.

Type:

list of arrays of complex

realdisp_associated_displacements

The lists of complex displacements that can be simultaneously computed from this set of real space displacement.

Type:

list of tuples

set_results(job_handler, data_type, prefix=None, fill=0.0, config_file='finite_displacements.yml')

Get finite displacements results from calculations.

Parameters:
  • job_handler (ComputeJobSeries or JobsDB) – Either a ComputeJobSeries or a JobsDB object as the interface to the jobs.

  • data_type (str, choices of "forces" and "energy") – Type of data to read from compute jobs.

  • prefix (str, optional, default to None) – Prefix to the jobnames.

  • config_file (str, optional, default to "finite_displacements.yml") – The name of the finite displacements configuration file.

_fp_results

The result of FrozenPhonons calculations.

Type:

list of arrays of complex, length of number of measurements

class principia_materia.phonon_id.frozen_phonons.LoneQ_FP(structure, Qpoint, supa=None, hidden_order=1, fdtype='c', tol=1e-06)

Bases: FrozenPhonons

LoneQ approach to comphe phonon and their interactions.

This approach does not consider symmetry, and simply find combinations of the reciprocal naive displacements basis and compute the phonon interactions for each unique dynamic tensor elements..

set_displacements()

Abstract method to construct displacements vectors for finite displacements.

set_displacements_basis()

Abstract method to construct the displacements basis vectors.

set_dynamic_tensor()

Construct the dynamic tensor after finite displacements result is computed.

_dynamic_tensor

The dynamic tensor at the Q-point.

Type:

array of complex, shape(norbitals, ) * order

set_naive_displacements()

Find the unique naive displacements to compute phonon interactions.

set_naive_displacements_basis()

Construct the naive reciprocal space displacement basis.

The displacement amplitudes follow the convention discussed in the paper Phys. Rev. B 100, 014303 (2019).

\[u_{\textbf{t}}*{\left(b, \beta\right)} = \sum_{\textbf{q} \in \hat{S}_{BZ}} u_{\textbf{q}}*{\left(b, \beta\right)} e^{2\pi\imag \textbf{t} \cdot \textbf{q}}\]
principia_materia.phonon_id.frozen_phonons.find_realspace_displacements(Qpoint, displabels, allreal=False, tol=1e-06)

Find realspace displacements from complex displacements basis.

In Lone-Q and Lone-ID finite displacements approaches, we find displacements in the reciprocal space first, these displacements are often complex numbers. However, we are only able to displacement atoms in real space. Thus, we need to find minumum the real/imag pairs of these displacements to be able to compute the derivatives.

Parameters:
  • Qpoint (array of Fraction, shape (order-hidden_order, dim)) – The Q-point.

  • displabels (tuple of tuples) – The reciprocal space displacement basis labels for a measurement.

  • allreal (bool, optional, default to False) – Whether the displacements are wll real or not.

  • tol (float) – Error tolerence.

Returns:

  • cs_pairs (array of int) – The cos/sin pairs (real/imag pairs) that are required to compute the complex displacements of the measurement. The array is in row convention, every row is a real space displacement denoted with a cos/sin pair. There are 2 possible integers in the array: 0 for cos(real) and 1 for sin(imag).

  • transformation_matrix (array of complex) – The transformation matrix from the complex displacement derivatives to the real space displacement derivatives. Each row of the matrix is the coefficients of the linear combination of the complex displacements that make the real space displacement.

  • associated_displacements (list of tuples) – A list of complex displacements that can be simultaneously computed from this set of real space displacement.

principia_materia.phonon_id.frozen_phonons.find_unique(input_list, return_index=False)

Find unique elements in a list.

Parameters:
  • input_list (list) – An input list.

  • return_index (bool, optional, default False) – Whether to return the indices of the unique elements.

principia_materia.phonon_id.frozen_phonons.format_irrep_instances(irrep_instances)
principia_materia.phonon_id.frozen_phonons.format_measurement_label(measurement_index, cs_indicies=None, zero_padding=0)

Format measurement indices into str labels.

Parameters:
  • measurement_index (int) – The index of the measurement.

  • cs_indicies (array of int, optional, default to None) – Array denoting the cos/sin pairs of the measurement.

  • zero_padding (int, optional, default to 0) – The amount of zero padding on the left of the measurement indices.

Returns:

label – The label for the measurement.

Return type:

str

principia_materia.phonon_id.frozen_phonons.get_FrozenPhonons_hdf5_wrapper()

HDF5 data wrapper for FrozenPhonons class.

principia_materia.phonon_id.frozen_phonons.get_LoneQ_from_hdf5(h5file, root_directory=None, fdseries_filename='fdseries.hdf5', tol=1e-06)

Read data from HDF5 file and construct a LoneQ_FP object.

Parameters:
  • h5file (str or h5py.File/h5py.Group obejct) – Path to a HDF5 file or a h5py.File/h5py.Group object.

  • tol (float, optional, default to 1.0E-6) – Error tolerance.

principia_materia.phonon_id.frozen_phonons.get_dynamic_tensor_wrapper()

HDF5 data wrapper for dynamic tensor of a q-point.

principia_materia.phonon_id.frozen_phonons.is_qpoint_real(qpoint)
principia_materia.phonon_id.frozen_phonons.load_fp_errortails(filename)

Load errortail data of frozen phonon runs from HDF5 file.

Parameters:

filename (str) – Path to HDF5 file.

principia_materia.phonon_id.frozen_phonons.make_displacement_identifier(Qpoint, displabels, cs=None, sort=True)

Make the nested tuple displacement identifier for a measurement.

Parameters:
  • Qpoint (array of Fractions) – A Q-point, an array of q-points.

  • displabels (list of tuples or str) – The labels for the displacement.

  • cs (list of str or int, optional, default to None) – The cos/sin labels for the displaement.

  • sort (bool, optional, default to True) – Whether to sort the q-point & displacement label tuple, when making the identifier.

Returns:

displacement_identifier – The identifier for the displacement.

Return type:

nested tuple

principia_materia.phonon_id.frozen_phonons.save_LoneQ_to_hdf5(obj, h5file='loneq_fp.hdf5', include_fdseires=False)

Save data of a LoneQ_FP object into a HDF5 file.

Parameters:
  • obj (LoneQ_FP) – A LoneQ_FP object.

  • h5file (str or h5py.File/h5py.Group obejct, optional, default to "loneq_fp.hdf5") – Path to a HDF5 file or a h5py.File/h5py.Group object.

principia_materia.phonon_id.frozen_phonons.save_dynamic_tensor_to_hdf5(obj, h5file='dynamic_tensor.hdf5')

Save data of a dynamic tensor of a q-poin into a HDF5 file.

Parameters:
  • obj (FrozenPhonons) – A FrozenPhonons object.

  • h5file (str or h5py.File/h5py.Group obejct, optional, default to "dynamic_tensor.hdf5") – Path to a HDF5 file or a h5py.File/h5py.Group object.

principia_materia.phonon_id.frozen_phonons.save_fp_errortails(output, fdtype, deltas, fd_values, result, xcoef, pick, penalty, overwrite=True)

Save errortail data of frozen phonon runs into HDF5 file.