Setting tags¶
Most of the setting tags have respective command-line options (Command options). When both of equivalent command-line option and setting tag are set simultaneously, the command-line option supersedes the setting tag. The configuration file is recommended to place at the first position for the mixed use of setting tags and command-line options, i.e.,
phonopy setting.conf [command-line-options]
For specifying real and reciprocal points, fractional values
(e.g. 1/3
) are accepted. However fractional values must not
have space among characters (e.g. 1 / 3
) are not allowed.
Basic tags¶
DIM
¶
The supercell is created from the input unit cell. When three integers are specified, a supercell elongated along axes of unit cell is created.
DIM = 2 2 3
In this case, a 2x2x3 supercell is created.
When nine integers are specified, the supercell is created by multiplying the supercell matrix \(M_\mathrm{s}\) with the unit cell. For example,
DIM = 0 1 1 1 0 1 1 1 0
the supercell matrix is
where the rows correspond to the first three, second three, and third three sets of numbers, respectively. When lattice parameters of unit cell are the column vectors of \(\mathbf{a}_\mathrm{u}\), \(\mathbf{b}_\mathrm{u}\), and \(\mathbf{c}_\mathrm{u}\), those of supercell, \(\mathbf{a}_\mathrm{s}\), \(\mathbf{b}_\mathrm{s}\), \(\mathbf{c}_\mathrm{s}\), are determined by,
Be careful that the axes in POSCAR
is defined by three row
vectors, i.e., \(( \mathbf{a}_\mathrm{u} \; \mathbf{b}_\mathrm{u}
\; \mathbf{c}_\mathrm{u} )^T\).
PRIMITIVE_AXES
or PRIMITIVE_AXIS
¶
When specified, transformation from the input unit cell to the primitive cell is performed. With this, the primitive cell basis vectors are used as the coordinate system for the phonon calculation. The transformation matrix is specified by nine values. The first, second, and third three values give the rows of the 3x3 matrix as follows:
PRIMITIVE_AXES = 0.0 0.5 0.5 0.5 0.0 0.5 0.5 0.5 0.0
Likewise,
PRIMITIVE_AXES = 0 1/2 1/2 1/2 0 1/2 1/2 1/2 0
The primitive cell for building the dynamical matrix is created by multiplying primitive-axis matrix \(\mathrm{M}_\mathrm{p}\). Let the matrix as,
where the rows correspond to the first three, second three, and third three sets of numbers, respectively.
When lattice parameters of unit cell (set by POSCAR
) are the
column vectors of \(\mathbf{a}_\mathrm{u}\),
\(\mathbf{b}_\mathrm{u}\), and \(\mathbf{c}_\mathrm{u}\),
those of supercell, \(\mathbf{a}_\mathrm{p}\),
\(\mathbf{b}_\mathrm{p}\), \(\mathbf{c}_\mathrm{p}\), are
determined by,
\(\mathrm{M}_\mathrm{p}\) is a change of basis matrix and so
\(\mathrm{M}_\mathrm{p}^{-1}\) must be an integer matrix. Be careful that
:math:the axes in POSCAR
is defined by three row vectors, i.e.,
\(( \mathbf{a}_\mathrm{u} \; \mathbf{b}_\mathrm{u} \;
\mathbf{c}_\mathrm{u} )^T\).
New in v1.14.0 PRIMITIVE_AXES = AUTO
is supported. This
enables to choose the transformation matrix automatically. Since the choice
of the primitive cell is arbitrary, it is recommended to use
PRIMITIVE_AXES = AUTO
to check if a possible transformation matrix
exists or not.
ATOM_NAME
¶
Chemical symbols
ATOM_NAME = Si O
The number of chemical symbols have to be same as that of the numbers
in the sixth line of POSCAR
.
Chemical symbols read by phonopy are overwritten by those written in
POSCAR
. See POSCAR
examples. In WIEN2k mode,
you don’t need to set this tag, i.e., chemical symbols are read from
the structure file.
EIGENVECTORS
¶
When this tag is ‘.TRUE.’, eigenvectors are calculated. With -p
option, partial density of states are calculated.
MASS
¶
This tag is not necessary to use usually, because atomic masses are automatically set from the chemical symbols.
Atomic masses of a primitive cell are overwritten by the values
specified. The order of atoms in the primitive cell that is defined by
PRIMITIVE_AXIS
tag can be shown using -v
option. It must be
noted that this tag does not affect to the symmetry search.
For example, when there are six atoms in a primitive cell, MASS
is
set as follows
MASS = 28.085 28.085 16.000 16.000 16.000 16.000
MAGMOM
¶
Symmetry of spin such as collinear magnetic moments is considered
using this tag. The number of values has to be equal to the number of
atoms in the unit cell, not the primitive cell or supercell. If this
tag is used with -d
option (CREATE_DISPLACEMENTS
tag),
MAGMOM
file is created. This file contains the MAGMOM
information of the supercell used for VASP. Unlike MAGMOM
in VASP,
*
can not be used, i.e., all the values (the same number of times
to the number of atoms in unit cell) have to be explicitly written.
MAGMOM = 1.0 1.0 -1.0 -1.0
FREQUENCY_CONVERSION_FACTOR
¶
Unit conversion factor of frequency from input values to your favorite unit can be specified, but the use should be limited and it is recommended to use this tag to convert the frequency unit to THz in some exceptional case, for example, a special force calculator whose physical unit system is different from the default setting of phonopy is used. If the frequency unit is different from THz, though it works just for seeing results of frequencies, e.g., band structure or DOS, it doesn’t work for derived values like thermal properties and mean square displacements.
The default values for calculators are those to convert frequency units to THz. The default conversion factors are shown at Default unit conversion factor of phonon frequency to THz. These are determined following the physical unit systems of the calculators. How to calcualte these conversion factors is explained at Physical unit conversion.
Displacement creation tags¶
CREATE_DISPLACEMENTS
¶
Supercells with displacements are created. This tag is used as the post process of phonon calculation.
CREATE_DISPLACEMENTS = .TRUE.
DIM = 2 2 2
DISPLACEMENT_DISTANCE
¶
Finite atomic displacement distance is set as specified value when
creating supercells with displacements. The default displacement
amplitude is 0.01 Angstrom, but when the wien2k
, abinit
,
Fleur
or turbomole
option is specified, the default value
is 0.02 Bohr.
DIAG
¶
When this tag is set .FALSE.
, displacements in diagonal directions
are not searched, i.e. all the displacements are along the lattice
vectors. DIAG = .FALSE.
is recommended if one of the lattice
parameter of your supercell is much longer or much shorter than the
other lattice parameters.
PM
¶
This tag specified how displacements are found. When PM = .FALSE.
,
least displacements that can calculate force constants are found. This
may cause less accurate result. When PM = .TRUE.
, all the
displacements that are opposite directions to the least displacements
are also found, which is called plus-minus displacements here. The
default setting is PM = AUTO
. Plus-minus displacements are
considered with this tag. If the plus and minus displacements are
symmetrically equivalent, only the plus displacement is found. This
may be in between .FALSE.
and .TRUE.
. You can check how it
works to see the file DISP
where displacement directions on atoms
are written.
RANDOM_DISPLACEMENTS
¶
The number of random displacement supercells are created by the specified positive integer values. In each supercell, all atoms are displaced in random direction with a constant displacement distance specified by DISPLACEMENT_DISTANCE tag. The random seed can be specified by RANDOM_SEED tag.
To obtain force constants with random displacements and respective forces, external force constants calculator is necessary.
CREATE_DISPLACEMENTS = .TRUE.
DIM = 2 2 2
RANDOM_DISPLACEMENTS = 20
DISPLACEMENT_DISTANCE = 0.03
RANDOM_SEED
¶
The random seed used for creating random displacements by RANDOM_DISPLACEMENTS tag. The value has to be 32bit unsigned int. The random seed is useful for crating the same random displacements with using the same number.
Band structure tags¶
BAND
and BAND_POINTS
¶
BAND
gives sampling band paths. The reciprocal points are
specified in reduced coordinates. The given points are connected for
defining band paths. When comma ,
is inserted between the points,
the paths are disconnected.
BAND_POINTS
gives the number of sampling points including the path
ends. The default value is BAND_POINTS = 51
.
An example of three paths, (0,0,0) to (1/2,0,1/2), (1/2,1/2,1) to (0,0,0), and (0,0,0) to (1/2,1/2,1/2), with 101 sampling points of each path are as follows:
BAND = 0 0 0 1/2 0 1/2, 1/2 1/2 1 0 0 0 1/2 1/2 1/2
BAND_POINTS = 101
BAND_LABELS
¶
Labels specified are depicted in band structure plot at the points of
band segments. The number of labels has to correspond to the number of
band paths specified by BAND
plus one. When LaTeX math style
expression such as \(\Gamma\) (\Gamma
) is expected, it is
probably necessary to place it between two $ characters.
BAND = 1/2 0 1/2 0 0 0 1/2 1/2 1/2
BAND_LABELS = X $\Gamma$ L
The colors of curves are automatically determined by matplotlib. The same color in a band segment shows the same kind of band. Between different band segments, the correspondence of colors doesn’t mean anything.
BAND_CONNECTION
¶
With this option, band connections are estimated from eigenvectors and
band structure is drawn considering band crossings. In sensitive
cases, to obtain better band connections, it requires to increase
number of points calculated in band segments by the BAND_POINTS
tag.
BAND = 1/2 0 1/2 0 0 0 1/2 1/2 1/2
BAND_POINTS = 101
BAND_CONNECTION = .TRUE.
Mesh sampling tags¶
Mesh sampling tags are used commonly for calculations of thermal properties and density of states.
MESH
, MP
, or MESH_NUMBERS
¶
MESH
numbers give uniform meshes in each axis. As the default
behavior, the center of mesh is determined by the Monkhorst-Pack
scheme, i.e., for odd number, a point comes to the center, and for
even number, the center is shifted half in the distance between
neighboring mesh points.
Examples of an even mesh with \(\Gamma\) center in two ways,
MESH = 8 8 8
GAMMA_CENTER = .TRUE.
MESH = 8 8 8
MP_SHIFT = 1/2 1/2 1/2
If only one float value is given, e.g., MESH = 100.0
,
\(\Gamma\) centred sampling mesh is generated with the mesh
numbers \((N_{\mathbf{a}^*}, N_{\mathbf{b}^*}, N_{\mathbf{c}^*})\)
computed following the convention of the VASP automatic
k-point generation, which is
where \(l\) is the value to be specified. With this,
GAMMA_CENTER
becomes simply ignored, but MP_SHIFT
works on top
of the \(\Gamma\) centred sampling mesh.
MESh = 100
MP_SHIFT
¶
MP_SHIFT
gives the shifts in direction along the corresponding
reciprocal axes (\(a^*\), \(b^*\), \(c^*\)). 0 or 1/2
(0.5) can be used as these values. 1/2 means the half mesh shift with
respect to neighboring grid points in each direction.
GAMMA_CENTER
¶
Instead of employing the Monkhorst-Pack scheme for the mesh sampling,
\(\Gamma\) center mesh is used. The default value is .FALSE.
.
GAMMA_CENTER = .TRUE.
WRITE_MESH
¶
With a dense mesh, with eigenvectors, without mesh symmetry, sometimes
its output file mesh.yaml
or mesh.hdf5
can be huge. However
when those files are not needed, e.g., in (P)DOS calculation,
WRITE_MESH = .FALSE.
can disable to write out those files. With
(P)DOS calculation, DOS output files are obtained even with
WRITE_MESH = .FALSE.
. The default setting is .TRUE.
.
WRITE_MESH = .FALSE.
Phonon density of states (DOS) tags¶
Phonon density of states (DOS) is calcualted either with a linear tetrahedron method (default) or smearing method. Phonons are calculated on a sampling mesh, therefore these tags must be used with Mesh sampling tags. The physical unit of horizontal axis is that of frequency that the user employs, e.g., THz, and that of vertical axis is {no. of states}/({unit cell} x {unit of the horizontal axis}). If the DOS is integrated over the frequency range, it will be \(3N_\mathrm{a}\) states, where \(N_\mathrm{a}\) is the number of atoms in the unit cell.
Phonon-DOS is formally defined as
where \(N\) is the number of unit cells and \(\lambda = (\nu, \mathbf{q})\) with \(\nu\) as the band index and \(\mathbf{q}\) as the q-point. This is computed on a set of descritized sampling frequency points for which \(\omega\) is specified arbitrary using DOS_RANGE. The phonon frequencies \(\omega_\lambda\) are obtained on a sampling mesh whose the number of grid points being \(N\). In the smearing method, the delta function is replaced by normal distribution (Gaussian function) with the standard deviation specified by SIGMA. In the tetrahedron method, the Brillouin integration is made analytically within tetrahedra in reciprocal space.
DOS
¶
This tag enables to calculate DOS. This tag is automatically set when
PDOS
tag or -p
option.
DOS = .TRUE.
DOS_RANGE
¶
DOS_RANGE = 0 40 0.1
Total and partial density of states are drawn with some parameters. The example makes DOS be calculated from frequency=0 to 40 with 0.1 pitch.
FMIN, FMAX, and FPITCH can be alternatively used to specify the minimum and maximum frequencies (the first and second values).
FMIN
, FMAX
, and FPITCH
¶
The uniform frequency sampling points for phonon-DOS calculation are
specified. FMIN
and FMAX
give the minimum, maximum frequencies
of the range, respectively, and FPITCH
gives the frequency pitch
to be sampled. These three values are the same as those that can be
specified by DOS_RANGE
.
PDOS
¶
Projected DOS is calculated using this tag. The formal definition is written as
where \(j\) is the atom indices and \(\hat{\mathbf{n}}\) is the unit projection direction vector. Without specifying PROJECTION_DIRECTION or XYZ_PROJECTION, PDOS is computed as sum of \(g^j(\omega, \hat{\mathbf{n}})\) projected onto Cartesian axes \(x,y,z\), i.e.,
The atom indices \(j\) are specified by
PDOS = 1 2, 3 4 5 6
These numbers are those in the primitive cell. ,
separates the
atom sets. In this example, atom 1 and 2 are summarized as one curve
and atom 3, 4, 5, and, 6 are summarized as another curve.
PDOS = AUTO
is supported To group symmetrically equivalent atoms
automatically.
PDOS = AUTO
EIGENVECTORS = .TRUE.
and MESH_SYMMETRY = .FALSE.
are
automatically set, therefore the calculation takes much more time than
usual DOS calculation. With a very dense sampling mesh, writing data
into mesh.yaml
or mesh.hdf5
can be unexpectedly huge. If only
PDOS is necessary but these output files are unnecessary, then it is
good to consider using WRITE_MESH = .FALSE.
(WRITE_MESH).
PROJECTION_DIRECTION
¶
Eigenvectors are projected along the direction specified by this tag. Projection direction is specified in reduced coordinates, i.e., with respect to a, b, c axes.
PDOS = 1, 2
PROJECTION_DIRECTION = -1 1 1
XYZ_PROJECTION
¶
PDOS is calculated using eigenvectors projected along x, y, and z
Cartesian coordinates. The format of output file projected_dos.dat
becomes different when using this tag, where phonon-mode-frequency and
x, y, and z components of PDOS are written out in the order:
frequency atom1_x atom1_y atom1_z atom2_x atom2_y atom2_z ...
With -p
option, three curves are drawn. These correspond to
sums of all projections to x, sums of all projections to y, and sums
of all projections to z composents of eigenvectors, respectively.
XYZ_PROJECTION = .TRUE.
SIGMA
¶
A smearing method is used instead of a linear tetrahedron method. This tag also specifies the smearing width. The unit is same as that used for phonon frequency. The default value is the value given by the difference of maximum and minimum frequencies divided by 100.
SIGMA = 0.1
DEBYE_MODEL
¶
By setting .TRUE.
, DOS at lower phonon frequencies are fit to a
Debye model. By default, the DOS from 0 to 1/4 of the maximum phonon
frequencies are used for the fitting. The function used to the fitting
is \(D(\omega)=a\omega^2\) where \(a\) is the parameter and
the Debye frequency is \((9N/a)^{1/3}\) where \(N\) is the
number of atoms in unit cell. Users have to unserstand that this is
not a unique way to determine Debye frequency. Debye frequency is
dependent on how to parameterize it.
DEBYE_MODEL = .TRUE.
MOMEMT
and MOMENT_ORDER
¶
Phonon moments for DOS and PDOS defined below are calculated using
these tags up to arbitrary order. The order is specified with
MOMENT_ORDER
(\(n\) in the formula). Unless MOMENT_ORDER
specified, the first and second moments are calculated.
The moments for DOS are given as
The moments for PDOS are given as
\(\omega_\text{min}\) and \(\omega_\text{max}\) are specified :using ref:dos_fmin_fmax_tags tags. When these are not specified, the moments are computed with the range of \(\epsilon < \omega < \infty\), where \(\epsilon\) is a small positive value. Imaginary frequencies are treated as negative real values in this computation, therefore it is not a good idea to set negative \(\omega_\text{min}\).
MOMENT = .TRUE.
MOMENT_ORDER = 3
Thermal displacements¶
TDISP
, TMAX
, TMIN
, and TSTEP
¶
Mean square displacements projected to Cartesian axes as a function of
temperature are calculated from the number of phonon excitations. The
usages of TMAX
, TMIN
, TSTEP
tags are same as those in
thermal properties tags. Phonon
frequencies in THz, which is the default setting of phonopy, are used to
obtain the mean square displacements, therefore physical units have to
be set properly for it (see Interfaces to calculators.) The result
is given in \(\text{Angstrom}^2\) and writen into
thermal_displacements.yaml
. See the detail of the method,
Thermal displacement. These tags must be used with
Mesh sampling tags
Optionally, FMIN
tag (--fmin
option) with a small value is
recommened to be set when q-points at \(\Gamma\) point or near
\(\Gamma\) point (e.g. using very dense sampling mesh) are sampled
to avoid divergence. FMAX
tag (--fmax
option) can be used to
specify an upper bound of phonon frequencies where the phonons are
considered in the summation. The projection is applied along arbitrary
direction using PROJECTION_DIRECTION
tag
(PROJECTION_DIRECTION).
mesh.yaml
or mesh.hdf5
is not written out from phonopy-1.11.14.
TDISP = .TRUE.
PROJECTION_DIRECTION = 1 1 0
TDISPMAT
, TMAX
, TMIN
, and TSTEP
¶
Mean square displacement matricies are calculated. The definition is
shown at Thermal displacement. Phonon frequencies in THz, which
is the default setting of phonopy, are
used to obtain the mean square displacement matricies, therefore
physical units have to be set properly for it (see
Interfaces to calculators.) The result is given in
\(\text{Angstrom}^2\) and writen into
thermal_displacement_matrices.yaml
where six matrix elements are
given in the order of xx, yy, zz, yz, xz, xy. In this yaml file,
displacement_matrices
and displacement_matrices_cif
correspond
to \(\mathrm{U}_\text{cart}\) and \(\mathrm{U}_\text{cif}\)
defined at Mean square displacement matrix, respectively.
Optionally, FMIN
tag (--fmin
option) with a small value is
recommened to be set when q-points at \(\Gamma\) point or near
\(\Gamma\) point (e.g. using very dense sampling mesh) are sampled
to avoid divergence. FMAX
tag (--fmax
option) can be used to
specify an upper bound of phonon frequencies where the phonons are
considered in the summation.
The 3x3 matrix restricts distribution of each atom around the equilibrium position to be ellipsoid. But the distribution is not necessarily to be so.
mesh.yaml
or mesh.hdf5
is not written out from phonopy-1.11.14.
TDISPMAT = .TRUE.
TDISPMAT_CIF
¶
This tag specifis a temperature (K) at which thermal displacement is
calculated and the mean square displacement matrix is written to the
cif file tdispmat.cif
with the dictionary item aniso_U
. Phonon
frequencies in THz, which is the default setting of phonopy, are used
to obtain the mean square displacement matricies, therefore physical
units have to be set properly for it (see
Interfaces to calculators.) The result is given in
\(\textrm{Angstrom}^2\).
mesh.yaml
or mesh.hdf5
is not written out from phonopy-1.11.14.
TDISPMAT_CIF = 1273.0
Specific q-points¶
QPOINTS
¶
When q-points are supplied, those phonons are calculated. Q-points are specified successive values separated by spaces and collected by every three values as vectors in reciprocal reduced coordinates.
QPOINTS = 0 0 0 1/2 1/2 1/2 1/2 0 1/2
With QPOINTS = .TRUE.
, q-points are read from QPOITNS
file
(see the file format at QPOINTS) in curret directory
phonons at the q-points are calculated.
QPOINTS = .TRUE.
WRITEDM
¶
WRITEDM = .TRUE.
Dynamical matrices \(D\) are written into qpoints.yaml
in the following \(6N\times3N\) format, where N is the number of atoms in
the primitive cell.
The physical unit of dynamical matrix is [unit of force] / ([unit of
displacement] * [unit of mass])
, i.e., square of the unit of phonon
frequency before multiplying the unit conversion factor
(see FREQUENCY_CONVERSION_FACTOR).
and \(D_{jj'}\) is
where j and j’ are the atomic indices in the primitive cell. The
phonon frequencies may be recovered from qpoints.yaml
by writing a
simple python script. For example, qpoints.yaml
is obtained for
NaCl at \(q=(0, 0.5, 0.5)\) by
phonopy --dim="2 2 2" --pa="0 1/2 1/2 1/2 0 1/2 1/2 1/2 0" --qpoints="0 1/2 1/2" --writedm
and the dynamical matrix may be used as
#!/usr/bin/env python
import yaml
import numpy as np
data = yaml.load(open("qpoints.yaml"))
dynmat = []
dynmat_data = data['phonon'][0]['dynamical_matrix']
for row in dynmat_data:
vals = np.reshape(row, (-1, 2))
dynmat.append(vals[:, 0] + vals[:, 1] * 1j)
dynmat = np.array(dynmat)
eigvals, eigvecs, = np.linalg.eigh(dynmat)
frequencies = np.sqrt(np.abs(eigvals.real)) * np.sign(eigvals.real)
conversion_factor_to_THz = 15.633302
print frequencies * conversion_factor_to_THz
Non-analytical term correction¶
NAC
¶
Non-analytical term correction is applied to dynamical
matrix. BORN
file has to be prepared in the current directory. See
BORN (optional) and Non-analytical term correction.
The default method is NAC_METHOD = GONZE
after v1.13.0.
NAC = .TRUE.
NAC_METHOD
¶
The method of non-analytical term correction is chosen by this tag
between two, NAC_METHOD = GONZE
(Correction by dipole-dipole interaction) and
NAC_METHOD = WANG
(Interpolation scheme at general q-points with non-analytical term correction), and the default is
the former after v1.13.0.
Q_DIRECTION
¶
This tag is used to activate non-analytical term correction (NAC) at
\(\mathbf{q}\rightarrow\mathbf{0}\), i.e. practically
\(\Gamma\)-point, because NAC is direction
dependent. With this tag, \(\mathbf{q}\) is specified in the
fractional coordinates of the reciprocal basis vectors. Only the
direction has the meaning. Therefore Q_DIRECTION = 1 1 1
and
Q_DIRECTION = 2 2 2
give the same result. This tag is valid for
QPOINTS
, IRREPS
, and MODULATION
tags.
Away from \(\Gamma\)-point, this setting is ignored and the specified q-point is used as the q-direction.
QPOINTS = 0 0 0
NAC = .TRUE.
Q_DIRECTION = 1 0 0
Group velocity¶
GROUP_VELOCITY
¶
Group velocities at q-points are calculated by using this tag. The group velocities are written into a yaml file corresponding to the run mode in Cartesian coordinates. The physical unit depends on physical units of input files and frequency conversion factor, but if VASP and the default settings (e.g., THz for phonon frequency) are simply used, then the physical unit will be Angstrom THz.
GROUP_VELOCITY = .TRUE.
Technical details are shown at Group velocity.
GV_DELTA_Q
¶
The reciprocal distance used for finite difference method is
specified. The default value is 1e-5
for the method of non-analytical
term correction by Gonze et al.. In other case, unless this tag is
specified, analytical derivative is used instead of the finite
difference method.
GV_DELTA_Q = 0.01
Symmetry¶
SYMMETRY_TOLERANCE
¶
This is used to set geometric tolerance to find symmetry of crystal structure. The default value is 1e-5. In general, it is not a good idea to loosen the tolerance. It is recommended to symmetrize crystal structure before starting phonon calculation, e.g., using --symmetry option.
SYMMETRY_TOLERANCE = 1e-3
MESH_SYMMETRY
¶
Symmetry search on the reciprocal sampling mesh is disabled by setting
MESH_SYMMETRY = .FALSE.
. In some case such as hexagonal systems or
primitive cells of cubic systems having F and I-centrings, the results
with and without mesh symmetry give slightly different values for
those proprerties that can employ mesh symmetry. This happens when the
uniform sampling mesh made along basis vectors desn’t have the same
crystallographic point group as the crystal itself. This symmetry
breaking may be also seen by the fact that weight
written in
mesh.yaml
can be different from possible order of product group of
site-symmetry group and time revesal symmetry. Generally the
difference becomes smaller when increasing the sampling mesh numbers.
FC_SYMMETRY
¶
Changed at v1.12.3
Previously this tag required a number for the iteration. From version 1.12.3, the way of symmetrization for translation invariance is modified and this number became unnecessary.
This tag is used to symmetrize force constants by translational
symmetry and permutation symmetry with .TRUE.
or .FALSE.
.
FC_SYMMETRY = .TRUE.
From the translation invariance condition,
where i and j are the atom indices, and \(\alpha\) and \(\beta\) are the Catesian indices for atoms i and j, respectively. When this condition is broken, the sum gives non-zero value. This value is subtracted from the diagonal blocks. Force constants are symmetric in each pair as
Mind that the other symmetries of force constants, i.e., the
symmetry from crystal symmetry or rotational symmetry, are broken to
use FC_SYMMETRY
.
Force constants¶
FORCE_CONSTANTS
¶
FORCE_CONSTANTS = READ
There are three values to be set, which are READ
and WRITE
,
and .FALSE.
. The default is .FALSE.
. When FORCE_CONSTANTS =
READ
, force constants are read from FORCE_CONSTANTS
file. With
FORCE_CONSTANTS = WRITE
, force constants calculated from
FORCE_SETS
are written to FORCE_CONSTANTS
file.
The file format of FORCE_CONSTANTS
is shown
here.
FULL_FORCE_CONSTANTS
¶
FULL_FORCE_CONSTANTS = .TRUE.
is used to compute full supercell
constants matrix. The default setting is .FALSE.
. By .TRUE.
or
.FALSE.
, the array shape becomes (n_patom, n_satom, 3, 3)
or
(n_satom, n_satom, 3, 3)
, respectively. The detail is found at
FORCE_CONSTANTS and force_constants.hdf5.
READ_FORCE_CONSTANTS
¶
READ_FORCE_CONSTANTS = .TRUE.
is equivalent to FORCE_CONSTANTS =
READ
.
WRITE_FORCE_CONSTANTS
¶
WRITE_FORCE_CONSTANTS = .TRUE.
is equivalent to FORCE_CONSTANTS =
WRITE
.
FC_CALCULATOR
¶
External force constants calculator can be used using this
tag. Currently ALM
is supported. The phonopy’s default force
constants calculator is based on finite difference method, for which
atomic displacements are made systematically. The following is the
list of the force constants calculator currently possible to be
invoked from phonopy.
ALM
¶
New in v2.3 ALM (https://github.com/ttadano/ALM) is based on fitting approach and any displacements set of atoms in supercell can be handled. For example, random displacements generated by RANDOM_DISPLACEMENTS can be used to compute force constants. To use ALM, its python module has to be installed. The installation instruction is found here.
When ALM is used, please cite the paper: T. Tadano and S. Tsuneyuki, J. Phys. Soc. Jpn. 87, 041015 (2018).
FC_CALCULATOR = ALM
Create animation file¶
ANIME_TYPE
¶
ANIME_TYPE = JMOL
There are V_SIM
, ARC
, XYZ
, JMOL
, and POSCAR
settings. Those may be viewed by v_sim
, gdis
, jmol
(animation), jmol
(vibration), respectively. For POSCAR
, a set
of POSCAR
format structure files corresponding to respective
animation images are created such as APOSCAR-000
,
APOSCAR-001
,….
There are several parameters to be set in the ANIME
tag.
ANIME
¶
The format of ``ANIME`` tag was modified after ver. 0.9.3.3.
For v_sim¶
ANIME = 0.5 0.5 0
The values are the q-point to be calculated. An animation file of
anime.ascii
is generated.
For the other animation formats¶
Phonon is only calculated at \(\Gamma\) point. So q-point is not necessary to be set.
anime.arc
, anime.xyz
, anime.xyz_jmol
, or APOSCAR-*
are generated according to the ANIME_TYPE
setting.
ANIME = 4 5 20 0.5 0.5 0
The values are as follows from left:
Band index given by ascending order in phonon frequency.
Magnitude to be multiplied. In the harmonic phonon calculation, there is no amplitude information obtained directly. The relative amplitude among atoms in primitive cell can be obtained from eigenvectors with the constraint of the norm or the eigenvectors equals one, i.e., number of atoms in the primitive is large, the displacements become small. Therefore this has to be adjusted to make the animation good looking.
Number of images in one phonon period.
(4-6) Shift of atomic points in reduced coordinate in real space. These values can be omitted and the default values are
0 0 0
.
For anime.xyz_jmol
, the first and third values are not used,
however dummy values, e.g. 0, are required.
Create modulated structure¶
MODULATION
¶
The MODULATION
tag is used to create a crystal structure with
displacements along normal modes at q-point in the specified supercell
dimension.
Atomic displacement of the j-th atom is created from the real part of the eigenvectors with amplitudes and phase factors as
where \(A\) is the amplitude, \(\phi\) is the phase, \(N_\mathrm{a}\) is the number of atoms in the supercell specified in this tag and \(m_j\) is the mass of the j-th atom, \(\mathbf{q}\) is the q-point specified, \(\mathbf{r}_{jl}\) is the position of the j-th atom in the l-th unit cell, and \(\mathbf{e}_j\) is the j-th atom part of eigenvector. Convention of eigenvector or dynamical matrix employed in phonopy is shown in Dynamical matrix.
If several modes are specified as shown in the example above, they are
overlapped on the structure. The output filenames are
MPOSCAR...
. Each modulated structure of a normal mode is written
in MPOSCAR-<number>
where the numbers correspond to the order of
specified sets of modulations. MPOSCAR
is the structure where all
the modulations are summed. MPOSCAR-orig
is the structure without
containing modulation, but the dimension is the one that is specified.
Some information is written into modulation.yaml
.
Usage¶
The first three (nine) values correspond to supercell dimension
(supercell matrix) like the DIM tag. The following
values are used to describe how the atoms are modulated. Multiple sets
of modulations can be specified by separating by comma ,
. In each
set, the first three values give a Q-point in the reduced coordinates
in reciprocal space. Then the next three values are the band index
from the bottom with ascending order, amplitude, and phase factor in
degrees. The phase factor is optional. If it is not specified, 0 is
used.
Before multiplying user specified phase factor, the phase of
the modulation vector is adjusted as the largest absolute value,
\(\left|\mathbf{e}_j\right|/\sqrt{m_j}\), of element of
3N dimensional modulation vector to be real. The complex modulation
vector is shown in modulation.yaml
.
MODULATION = 3 3 1, 1/3 1/3 0 1 2, 1/3 1/3 0 2 3.5
MODULATION = 3 3 1, 1/3 1/3 0 1 2, 1/3 0 0 2 2
MODULATION = 3 3 1, 1/3 1/3 0 1 1 0, 1/3 1/3 0 1 1 90
MODULATION = -1 1 1 1 -1 1 1 1 -1, 1/2 1/2 0 1 2
Characters of irreducible representations¶
IRREPS
¶
Characters of irreducible representations (Irreps) of phonon modes are
shown. For this calculation, a primitive cell has to be used. If the
input unit cell is a non-primitive cell, it has to be transformed to a
primitive cell using PRIMITIVE_AXES
tag.
The first three values gives a q-point in reduced coordinates to be calculated. The degenerated modes are searched only by the closeness of frequencies. The frequency difference to be tolerated is specified by the fourth value in the frequency unit that the user specified.
IRREPS = 0 0 0 1e-3
Symbols of Irreps for the 32 point group types at the \(\Gamma\) point are shown but not at non-\(\Gamma\) point.
SHOW_IRREPS
¶
Irreducible representations are shown along with character table.
IRREPS = 1/3 1/3 0
SHOW_IRREPS = .TRUE.
LITTLE_COGROUP
¶
Show irreps of little co-group (point-group of wavevector) instead of little group.
IRREPS = 0 0 1/8
LITTLE_COGROUP = .TRUE.
Input/Output file control¶
FC_FORMAT
, READFC_FORMAT
, WRITEFC_FORMAT
¶
There are two file-formats to store force constants. Currently
text style (TEXT
) and hdf5 (HDF5
)
formats are supported. The default file format is the text
style. Reading and writing force constants are
invoked by FORCE_CONSTANTS tag. Using
these tags, the input/output formats are switched.
FC_FORMAT
affects to both input and output, e.g.:
FORCE_CONSTANTS = WRITE
FC_FORMAT = HDF5
READFC_FORMAT
and WRITEFC_FORMAT
can be used to control
input and output formats separately, i.e., the following setting to
convert force constants format is possible:
READ_FORCE_CONSTANTS = .TRUE.
WRITE_FORCE_CONSTANTS = .TRUE.
WRITEFC_FORMAT = HDF5
BAND_FORMAT
, MESH_FORMAT
, QPOINTS_FORMAT
¶
There are two file-formats to write the results of band structure,
mesh, and q-points calculations. Currently YAML (YAML
) and hdf5
(HDF5
) formats are supported. The default file format is the YAML
format. The file format is changed as follows:
BAND_FORMAT = HDF5
MESH_FORMAT = HDF5
QPOINTS_FORMAT = HDF5
HDF5
¶
The following output files are written in hdf5 format instead of their
original formats (in parenthesis) by HDF5 = .TRUE.
. In addition,
force_constants.hdf5
is read with this tag.
force_constants.hdf5
(FORCE_CONSTANTS
)mesh.hdf5
(mesh.yaml
)band.hdf5
(band.yaml
)qpoints.hdf5
(qpoints.yaml
)
HDF5 = .TRUE.
force_constants.hdf5
¶
With --hdf5
option and FORCE_CONSTANTS = WRITE
(--writefc
), force_constants.hdf5
is written.
With --hdf5
option and FORCE_CONSTANTS = READ
(--readfc
),
force_constants.hdf5
is read.
mesh.hdf5
¶
In the mesh sampling calculations (see Mesh sampling tags),
calculation results are written into mesh.hdf5
but not into
mesh.yaml
. Using this option may reduce the data output size and
thus writing time when mesh.yaml
is huge, e.g., eigenvectors are
written on a dense sampling mesh.
qpoints.hdf5
¶
In the specific q-points calculations (QPOINTS),
calculation results are written into qpoints.hdf5
but not into
qpoints.yaml
. With WRITEDM, dynamical matrices are also
stored in qpoints.hdf5
. Using this option may be useful with large
set of q-points with including eigenvector or dynamical matrix output.
band.hdf5
¶
In the band structure calculations (Band structure tags),
calculation results are written into band.hdf5
but not into
band.yaml
.
summary
¶
The following data may be optionally included in the summary yaml file
called phonopy_disp.yaml
/phonopy.yaml
in addition to other file
output settings. This happens at the end of the pre/post-process (after
running the phonopy
script):
force constants
force sets
dielectric constant
born effective charge
displacements
[all]
Including all relevant data in a single output file allows for a human readable convenient file format.
force constants
¶
The --include-fc
flag or setting INCLUDE_FC = .TRUE.
will cause
the force constants (if available) to be written as an entry in the
yaml summary file. The written force constants will reflect the
required/available format used during processing. So if --full-fc
is
set the entire matrix will be written.
force sets
¶
The --include-fs
flag or setting INCLUDE_FS = .TRUE.
will cause
the force sets (if available) to be written as an entry in the yaml summary
file.
dielectric constant
and born effective charge
¶
The --include-born
flag or setting INCLUDE_BORN = .TRUE.
will cause
the born effective charges and dielectric tensor (if available) to be
written as an entry in the yaml summary file. The values will only be written
if non-analytical term correction is set with the --nac
flag or by
setting NAC = .TRUE.
.
This is more convenient than keeping track of the BORN
file created by the user.
displacements
¶
The --include-disp
flag or setting INCLUDE_DISP = .TRUE.
will cause
displacements data (if available) to be written as an entry in the yaml summary file.
This is set by default when the phonopy
script is run in displacements
mode.
all
¶
All available data covered by the other include
flags can be written to the yaml
summary file using the --include-all
flag or by setting INCLUDE_ALL = .TRUE.
.
Force constants are not stored when force sets are stored.