Additional functions

calorine.tools.analyze_data(data, max_lag=None)[source]

Carries out an extensive analysis of the data series.

Parameters:
  • data (ndarray) – data series to compute autocorrelation function for

  • max_lag (Optional[int]) – maximum lag between two data points, used for computing autocorrelation

Returns:

calculated properties of the data including, mean, standard deviation, correlation length and a 95% error estimate.

Return type:

dict

calorine.tools.get_autocorrelation_function(data, max_lag=None)[source]

Returns autocorrelation function.

The autocorrelation function is computed using pandas.Series.autocorr.

Parameters:
  • data (ndarray) – data series to compute autocorrelation function for

  • max_lag (Optional[int]) – maximum lag between two data points

Return type:

ndarray

Returns:

calculated autocorrelation function

calorine.tools.get_correlation_length(data)[source]

Returns estimate of the correlation length of data.

The correlation length is taken as the first point where the autocorrelation functions is less than \(\exp(-2)\). If the correlation function never drops below \(\exp(-2)\) np.nan is returned.

If the correlation length cannot be computed since the auto-correlation function is unconverged the function returns None.

Parameters:

data (ndarray) – data series for which to the compute autocorrelation function

Return type:

Optional[int]

Returns:

correlation length

calorine.tools.get_elastic_stiffness_tensor(structure, clamped=False, epsilon=0.001, **kwargs)[source]

Calculate and return the elastic stiffness tensor in units of GPa for the given structure in Voigt form.

Parameters:
  • structure (Atoms) – input structure; should be fully relaxed

  • clamped (bool) – if False (default) return the relaxed elastic stiffness tensor; if True return the clamped ion elastic stiffness tensor

  • epsilon (float) – magnitude of the applied strain

  • kwargs – keyword arguments forwarded to the relax_structure function used for relaxing the structure when computing the relaxed stiffness tensor; it should not be necessary to change the default for the vast majority of use cases; use with care

Return type:

Stiffness tensor in units of GPa

calorine.tools.get_error_estimate(data, confidence=0.95)[source]

Returns estimate of standard error \(\mathrm{error}\) with confidence interval.

\[\mathrm{error} = t_\mathrm{factor} * \mathrm{std}(\mathrm{data}) / \sqrt{N_s}\]

where \(t_{factor}\) is the factor corresponding to the confidence interval and \(N_s\) is the number of independent measurements (with correlation taken into account).

If the correlation length cannot be computed since the auto-correlation function is unconverged the function returns None.

Parameters:

data (ndarray) – data series for which to estimate the error

Return type:

Optional[float]

Returns:

error estimate

calorine.tools.get_force_constants(structure, calculator, supercell_matrix, kwargs_phonopy={}, kwargs_generate_displacements={})[source]

Calculates the force constants for a given structure using phonopy, which needs to be cited if this function is used for generating data for publication. The function returns a Phonopy object that can be used to calculate, e.g., the phonon dispersion, the phonon density of states as well as related quantities such as the thermal displacements and the free energy.

Parameters:
  • structure (Atoms) – structure for which to compute the phonon dispersion; usually this is a primitive cell

  • calculator (SinglePointCalculator) – ASE calculator to use for the calculation of forces

  • supercell_matrix (ndarray) – specification of supercell size handed over to phonopy; should be a tuple of three values or a matrix

  • kwargs_phonopy (Dict[str, Any]) – Expert option: keyword arguments used when initializing the Phonopy object; this includes, e.g., the tolerance used when determining the symmetry (symprec) and parameters for the non-analytical corrections (nac_params)

  • kwargs_generate_displacements (Dict[str, Any]) – Expert option: keyword arguments to be handed over to the generate_displacements method; this includes in particular the distance keyword, which specifies the magnitude of the atomic displacement imposed when calculating the force constant matrix

Return type:

Phonopy

calorine.tools.get_primitive_structure(structure, no_idealize=True, to_primitive=True, symprec=1e-05)[source]

Returns the primitive structure using spglib. This is a convenience interface to the standardize_cell() function of spglib that works directly with ase Atoms objects.

Parameters:
  • structure (Atoms) – Input atomic structure.

  • no_idealize (bool) – If True lengths and angles are not idealized.

  • to_primitive (bool) – If True convert to primitive structure.

  • symprec (float) – Tolerance imposed when analyzing the symmetry.

Return type:

Atoms

calorine.tools.get_rtc_from_hac(hac, V, T, dt)[source]

Returns the running thermal conductivity (RTC) in W/m/K using a heat auto-correlation (HAC) as input.

Parameters:
  • hac (ndarray) – The HAC \(\langle j(t)j(0)\rangle\) in units of eV3/amu as given by GPUMD

  • V (float) – Volume of cell in Å3

  • T (float) – Temperature in Kelvin

  • dt (float) – Time step in the HAC in ps running thermal conductivity in W/m/K

Return type:

ndarray

calorine.tools.get_spacegroup(structure, symprec=1e-05, angle_tolerance=-1.0, style='international')[source]

Returns the space group of a structure using spglib. This is a convenience interface to the get_spacegroup() function of spglib that works directly with ase Atoms objects.

Parameters:
  • structure (Atoms) – Input atomic structure.

  • symprec (float) – Tolerance imposed when analyzing the symmetry.

  • angle_tolerance (float) – Tolerance imposed when analyzing angles.

  • style (str) – Space group notation to be used. Can be 'international' for the interational tables of crystallography (ITC) style (Hermann-Mauguin and ITC number) or 'Schoenflies' for the Schoenflies notation.

Return type:

str

calorine.tools.get_wyckoff_sites(structure, map_occupations=None, symprec=1e-05, include_representative_atom_index=False)[source]

Returns the Wyckoff symbols of the input structure. The Wyckoff labels can be conveniently attached as an array to the structure object as demonstrated in the examples section below.

By default the occupation of the sites is part of the symmetry analysis. If a chemically disordered structure is provided this will usually reduce the symmetry substantially. If one is interested in the symmetry of the underlying structure one can control how occupations are handled. To this end, one can provide the map_occupations keyword argument. The latter must be a list, each entry of which is a list of species that should be treated as indistinguishable. As a shortcut, if all species should be treated as indistinguishable one can provide an empty list. Examples that illustrate the usage of the keyword are given below.

Parameters:
  • structure (Atoms) – Input structure. Note that the occupation of the sites is included in the symmetry analysis.

  • map_occupations (Optional[List[List[str]]]) – Each sublist in this list specifies a group of chemical species that shall be treated as indistinguishable for the purpose of the symmetry analysis.

  • symprec (float) – Tolerance imposed when analyzing the symmetry using spglib.

  • include_representative_atom_index (bool) – If True the index of the first atom in the structure that is representative of the Wyckoff site is included in the symbol. This is in particular useful in cases when there are multiple Wyckoff sites sites with the same Wyckoff letter.

Return type:

List[str]

Examples

Wyckoff sites of a hexagonal-close packed structure:

>>> from ase.build import bulk
>>> structure = bulk('Ti')
>>> wyckoff_sites = get_wyckoff_sites(structure)
>>> print(wyckoff_sites)
['2d', '2d']

The Wyckoff labels can also be attached as an array to the structure, in which case the information is also included when storing the Atoms object:

>>> from ase.io import write
>>> structure.new_array('wyckoff_sites', wyckoff_sites, str)
>>> write('structure.xyz', structure)

The function can also be applied to supercells:

>>> structure = bulk('GaAs', crystalstructure='zincblende', a=3.0).repeat(2)
>>> wyckoff_sites = get_wyckoff_sites(structure)
>>> print(wyckoff_sites)
['4a', '4c', '4a', '4c', '4a', '4c', '4a', '4c',
 '4a', '4c', '4a', '4c', '4a', '4c', '4a', '4c']

Now assume that one is given a supercell of a (Ga,Al)As alloy. Applying the function directly yields much lower symmetry since the symmetry of the original structure is broken:

>>> structure.set_chemical_symbols(
...        ['Ga', 'As', 'Al', 'As', 'Ga', 'As', 'Al', 'As',
...         'Ga', 'As', 'Ga', 'As', 'Al', 'As', 'Ga', 'As'])
>>> print(get_wyckoff_sites(structure))
['8g', '8i', '4e', '8i', '8g', '8i', '2c', '8i',
 '2d', '8i', '8g', '8i', '4e', '8i', '8g', '8i']

Since Ga and Al occupy the same sublattice, they should, however, be treated as indistinguishable for the purpose of the symmetry analysis, which can be achieved via the map_occupations keyword:

>>> print(get_wyckoff_sites(structure, map_occupations=[['Ga', 'Al'], ['As']]))
['4a', '4c', '4a', '4c', '4a', '4c', '4a', '4c',
 '4a', '4c', '4a', '4c', '4a', '4c', '4a', '4c']

If occupations are to ignored entirely, one can simply provide an empty list. In the present case, this turns the zincblende lattice into a diamond lattice, on which case there is only one Wyckoff site:

>>> print(get_wyckoff_sites(structure, map_occupations=[]))
['8a', '8a', '8a', '8a', '8a', '8a', '8a', '8a',
 '8a', '8a', '8a', '8a', '8a', '8a', '8a', '8a']
calorine.tools.relax_structure(structure, fmax=0.001, steps=500, minimizer='bfgs', constant_cell=False, constant_volume=False, scalar_pressure=0.0, **kwargs)[source]

Relaxes the given structure.

Parameters:
  • structure (Atoms) – Atomic configuration to relax.

  • fmax (float) – Stop relaxation if the absolute force for all atoms falls below this value.

  • steps (int) – Maximum number of relaxation steps the minimizer is allowed to take.

  • minimizer (str) – Minimizer to use; possible values: ‘bfgs’, ‘lbfgs’, ‘fire’, ‘gpmin’, ‘bfgs-scipy’.

  • constant_cell (bool) – If True do not relax the cell or the volume.

  • constant_volume (bool) – If True relax the cell shape but keep the volume constant.

  • kwargs – Keyword arguments to be handed over to the minimizer; possible arguments can be found in the ASE documentation https://wiki.fysik.dtu.dk/ase/ase/filters.html#the-frechetcellfilter-class.

  • scalar_pressure (float) – External pressure in GPa.

Return type:

None