\(\newcommand{\AA}{\text{Å}}\)
CRYSTALClear package
CRYSTALClear.adsorb
- sub_ads_indices(structure)
Identify the indices of substrate and adsorbate atoms in the given structure.
- Parameters:
structure (pymatgen.core.structure.Structure) – The input structure.
- Returns:
dict – A dictionary containing the indices of adsorbate and substrate atoms. The dictionary has two keys: - ‘adsorbate’: A list of indices corresponding to the adsorbate atoms. - ‘substrate’: A list of indices corresponding to the substrate atoms.
CRYSTALClear.calculate
- cry_ads_energy(e_full_system, e_substrate, e_adsorbate)
Calculate the adsorption energy of a system.
- Parameters:
e_full_system (float) – Total energy of the full system.
e_substrate (float) – Energy of the substrate.
e_adsorbate (float) – Energy of the adsorbate.
- Returns:
float – Adsorption energy calculated as the difference between the total energy of the full system and the sum of the energies of the substrate and adsorbate.
- cry_shrink(structure, spacing=0.2)
Determine the number of unit cells needed to achieve a desired spacing.
- Parameters:
structure (pymatgen.core.structure.Structure) – The input structure.
spacing (float) – The desired spacing between unit cells. Default is 0.2.
- Returns:
int – The number of unit cells (rounded up) needed to achieve the desired spacing.
CRYSTALClear.convert
Functions that do conversion between data / file formats
- cry_ase2gui(structure, gui_file=None, vacuum=None, symmetry=True)
Transform an ASE Structure object into a Pymatgen structure object and then a CRYSTAL structure (gui) object. Vacuum layer is set to 500 Angstrom as the default of CRYSTAL for low symmstry systems
- Parameters:
structure (ASE Structure) – ASE Structure object.
gui_file (str) – CRYSTAL gui / fort.34 file.
vacuum (float) – Vacuum distance. Unit: Angstrom. If none, set the length of non-periodic direction to 500 Angstrom. Low dimensional systems only.
symmetry (bool) – Perform symmetry analysis.
- cry_bands2pmg(band, output, labels=None)
Transform a CRYSTAL bands object into a Pymatgen bands object. No projection is available for now.
- Parameters:
band (str|Band) – CRYSTAL fort.25 or BAND.DAT file or CRYSTALClear Band object
output (str) – CRYSTAL properties output file
labels (list[str]) – K point labels to display in the band structure.
- Returns:
BandStructureSymmLine – Pymatgen band structure object.
- cry_gui2ase(gui, vacuum=None, **kwargs)
Transform a CRYSTAL structure (gui) file into an ASE atoms object.
- Parameters:
gui (str|Crystal_gui) – CRYSTAL gui / fort.34 file or CRYSTALClear gui object
vacuum (float) – Vacuum distance. Unit: Angstrom. If none, set the
pbc
attribute of ASE atoms object. Low dimensional systems only.**kwargs – Passed to ASE Atoms constructor
- Returns:
Atoms – ASE atoms object.
- cry_gui2cif(gui, cif_file_name, vacuum=None, **kwargs)
Read a CRYSTAL structure (gui) file and save a cif file. The CifWriter object of Pymatgen is called. By default,
symprec = 0.01
is used.- Parameters:
gui (str|Crystal_gui) – CRYSTAL gui / fort.34 file or CRYSTALClear gui object
cif_file_name (str) – Name (including path) of the cif file to be saved
vacuum (float) – Vacuum distance. Unit: Angstrom. If none, set the
pbc
attribute of Pymatgen atoms object. Low dimensional systems only.**kwargs – Passed to Pymatgen CifWriter.
- cry_gui2pmg(gui, vacuum=None, molecule=True)
Transform a CRYSTAL structure (gui) object into a Pymatgen Structure object.
- Parameters:
gui (str|Crystal_gui) – CRYSTAL gui / fort.34 file or CRYSTALClear gui object
vacuum (float) – Vacuum distance. Unit: Angstrom. If none, set the
pbc
attribute of Pymatgen object. Low dimensional systems only.molecule (bool) – Generate a Molecule Pymatgen object for 0D structures.
- Returns:
Structure or Molecule – Pymatgen Structure or Molecule object.
- cry_gui2xyz(gui, xyz_file_name, **kwargs)
Transform a CRYSTAL structure (gui) file into an XYZ file.
- Parameters:
xyz_file_name (str) – Name of the XYZ file to be saved.
gui (str) – CRYSTAL gui / fort.34 file.
**kwargs – Passed to Pymatgen XYZ object.
- cry_out2ase(output, vacuum=None, initial=False, **kwargs)
Transform a CRYSTAL output object into an ASE atoms object.
- Parameters:
output (str|Crystal_output) – Crystal output file or Crystal_output object
vacuum (float) – Vacuum distance. Unit: Angstrom. If none, set the
pbc
attribute of ASE atoms object. Low dimensional systems only.initial (bool) – Read the last geometry of the output file.
**kwargs – Passed to ASE Atoms constructor
- Returns:
Atoms – ASE atoms object.
- cry_out2cif(output, cif_file_name, vacuum=None, initial=False, **kwargs)
Read geometry from a CRYSTAL output file and save it as a CIF file. The CifWriter object of Pymatgen is called. By default,
symprec = 0.01
is used.- Parameters:
output (str|Crystal_output) – Crystal output file or Crystal_output object
cif_file_name (str) – Name (including path) of the CIF file to be saved.
vacuum (float) – Vacuum distance. Unit: Angstrom. If none, set the
pbc
attribute of Pymatgen atoms object. Low dimensional systems only.initial (bool) – Read the last geometry of the output file.
**kwargs – Passed to Pymatgen CifWriter.
- cry_out2pmg(output, vacuum=None, initial=False, molecule=True)
Read geometry from a CRYSTAL output file and save it as a pymatgen structure object.
- Parameters:
output (str|Crystal_output) – Crystal output file or Crystal_output object
vacuum (float) – Vacuum distance. Unit: Angstrom. If none, set the
pbc
attribute of Pymatgen object. Low dimensional systems only.initial (bool) – Read the last geometry of the output file.
molecule (bool) – Generate a Molecule Pymatgen object for 0D structures.
- Returns:
Structure – Pymatgen Structure object.
- cry_out2xyz(output, xyz_file_name, initial=False, **kwargs)
Read geometry from a CRYSTAL output file and save it as an XYZ file.
- Parameters:
output (str|Crystal_output) – Crystal output file or Crystal_output object
xyz_file_name (str) – Name (including path) of the XYZ file to be saved.
initial (bool) – Read the last geometry of the output file.
**kwargs – Passed to Pymatgen XYZ object.
- cry_pmg2gui(structure, gui_file=None, pbc=None, vacuum=None, symmetry=True, zconv=None, **kwargs)
Save a pymatgen Structure object into a CRYSTAL gui file. Vacuum layer is set to 500 Angstrom as the default of CRYSTAL for low symmstry systems.
- Parameters:
structure (Structure | Molecule) – Pymatgen Structure / Molecule object.
gui_file (str) – CRYSTAL gui / fort.34 file.
pbc (list) – 1*3 boolian list. Implements periodicity along x, y and z directions. If none, the code will read it from input structure.
vacuum (float) – Vacuum distance. Unit: Angstrom. If none, set the length of non-periodic direction to 500 Angstrom. Low dimensional systems only.
symmetry (bool) – Do symmetry analysis.
zconv (list[list[int, int]]) – 1st element: The index of atom; 2nd element: The new conventional atomic number.
**kwargs – Passed to Pymatgen SpacegroupAnalyzer object. Valid only if
symmetry=True
.
CRYSTALClear.crystal_io
Objects of input / output files of CRYSTAL. Methods to edit or substract data from corresponding files are provided.
- class Crystal_input
Bases:
Crystal_inputBASE
Crystal input object inherited from the Crystal_inputBASE object. For the basic set-ups of keywords, please refer to manuals there.
- geom_from_cif(file, zconv=None, keyword='EXTERNAL', pbc=[True, True, True], gui_name='fort.34', **kwargs)
Read geometry from cif file and put infomation to geom block, either as ‘EXTERNAL’ or ‘CRYSTAL’. CIF files with a single geometry only.
CIF with Symmetry is required.
Note
CIF files should have the ‘.cif’ extension.
- Parameters:
file (str) – CIF file.
keyword (str) – ‘EXTERNAL’ or ‘CRYSTAL’.
zconv (list[list[int, int]]) – 1st element: The index of atom; 2nd element: The new conventional atomic number. Atoms of the irreducible unit is required.
pbc (list[bool]) – Valid only if
keyword = EXTERNAL
. Force to remove periodic boundary conditions along x, y, z axis.gui_name (str) – Valid only if
keyword = EXTERNAL
. Gui file’s name.**kwargs – Passed to Pymatgen SpacegroupAnalyzer object.
- geom_from_pmg(struc, zconv=None, keyword='EXTERNAL', gui_name='fort.34', **kwargs)
Read geometry defined by PyMatGen structure object and put infomation into geom block, either as ‘EXTERNAL’ or ‘CRYSTAL’.
See
geom_from_cif
for definition of arguments.Note
Coordinates of corresponding atoms may not consistent with the original CIF file if ‘CRYSTAL’ is used, in which case coordinates of another symmetry equivalent atom is used.
- bs_from_bse(name, element, zconv=None)
Download basis set definitions from BSE. A wrapper of BasisSetBASE.from_bse.
- Parameters:
name (str) – Basis set’s name.
element (list[str] | list[int]) – List of elements.
zconv (list[int]) – If not none, use the conventional atomic number. Its length must be the same as element. Its sequence must be consistent with basis set’s
- bs_from_string(fmt='crystal')
Define basis set from a string. A wrapper of BasisSetBASE.from_string.
- Parameters:
bs_str (str)
fmt (str) – Format string. Consistent with BSE python API.
- bs_from_file(fmt='crystal')
Define a basis set from a file. A wrapper of BasisSetBASE.from_file.
- class Crystal_output(output=None)
Bases:
object
This class reads a CRYSTAL output and generates an object.
- read_cry_output(output_name)
Deprecated
- get_dielectric_tensor()
Extracts the dielectric tensor from the output.
- Returns:
list – Dielectric tensor values.
- get_eigenvectors()
Extracts eigenvectors from the output.
- get_dimensionality()
Gets the dimensionality of the system.
- Returns:
self.dimensionality (int) – Dimensionality of the system.
- get_symmops()
Return the symmetry operators
- Returns:
self.symmops (numpy.ndarray) – Symmetry operators
- get_geometry(initial=True, write_gui=False, gui_name=None, symmetry='pymatgen', **kwargs)
Get the geometry.
- Parameters:
initial (bool) – Read the initial or last gemetry. Useful in case of geometry optimization.
write_gui (bool) – If True, write .gui file
gui_name (str) – Valid only if
write_gui = True
. Gui file is named as ‘gui_name’. If None, use ‘basename.gui’. The basename is the same as output file.symmetry (str) – Valid only if
write_gui = True
. ‘pymatgen’ to use symmetry info from a pymatgen SpacegroupAnalyzer; ‘initial’ to use symmstry information on output file. If None, no symmstry. Otherwise it is taken from the existing gui file.**kwargs – Valid only if
write_gui = True
andsymmetry = 'pymatgen'
. Passed to Pymatgen SpacegroupAnalyzer object.
- Returns:
self.geometry (Structure | Molecule) – A pymatgen Structure or molecule object.
- get_last_geom(write_gui_file=True, symm_info='pymatgen')
Return the last optimised geometry.
- get_lattice(initial=True)
Returns the lattice of the system. Unit: Angstrom.
- Parameters:
initial (bool) – Read the initial or last lattice. Useful in case of geometry optimization.
- Returns:
self.lattice (np.ndarray) – Lattice of the system.
- get_reciprocal_lattice(initial=True)
Returns the reciprocal lattice of the system. Unit: Angstrom^-1. Reciprocal lattice is consistent with CRYSTAL: 2pi is added.
- Parameters:
initial (bool) – Read the initial or last lattice. Useful in case of geometry optimization.
- Returns:
self.reciprocal_lattice (np.ndarray) – Lattice of the system.
- property sg_number
CRYSTAL 0~3D space group number. Before geometry editing.
- property sg_symbol
CRYSTAL 0~3D 0~3D space group symbol. Before geometry editing.
- property n_atoms
Number of atoms. After geometry editing.
- property atom_symbols
Atom symbols. After geometry editing.
- property atom_numbers
Conventional atom numbers. After geometry editing.
- property atom_positions
Composite fractional / Cartesian atomic coordinates. Consistent with CRYSTAL definitions. 3D: Fractional; 2D: Frac, Frac, Cart; 1D Frac, Cart, Cart; 0D: Cart, Cart, Cart. After geometry editing (before optimization).
- property atom_positions_cart
Cartesian atomic coordinates. After geometry editing (before optimization).
- property atom_positions_frac
Fractional atomic coordinates. After geometry editing (before optimization).
- get_trans_matrix()
Get cell transformation matrix
- Returns:
self.trans_matrix (np.ndarray) – 3*3 array of supercell expansion matrix
- get_primitive_geometry(initial=True, write_gui=False, gui_name=None, symmetry='pymatgen', **kwargs)
Get the primitive geometry, reduced by cell transformation matrix inversed.
Note
This is not the standard ‘primitive cell’. This method returns geometry before supercell expansion keywords such as ‘SUPERCEL’ or ‘SCELPHONO’.
Conventional atomic numbers are not available.
- Parameters:
initial (bool) – Read the initial or last geometry. Useful in case of geometry optimization.
write_gui (bool) – If True, write .gui file
gui_name (str) – Valid only if
write_gui = True
. Gui file is named as ‘gui_name’. If None, use ‘basename.gui’. The basename is the same as output file.symmetry (str) – Valid only if
write_gui = True
. ‘pymatgen’ to use symmetry info from a pymatgen SpacegroupAnalyzer; ‘initial’ to use symmstry information on output file. If None, no symmstry. Otherwise it is taken from the existing gui file.**kwargs – Valid only if
write_gui = True
andsymmetry = 'pymatgen'
. Passed to Pymatgen SpacegroupAnalyzer object.
- Returns:
self.primitive_geometry (Structure | Molecule) – A pymatgen Structure or molecule object.
- get_primitive_lattice(initial=True)
Returns the primitive lattice of the system reduced by cell transformation matrix inverse. Unit: Angstrom.
- Parameters:
initial (bool) – Read the initial or last lattice. Useful in case of geometry optimization.
- Returns:
self.primitive_lattice (np.ndarray) – Lattice of the system.
- get_primitive_reciprocal_lattice(initial=False)
Returns the primitive reciprocal lattice of the system before expansion by cell transformation matrix inverse. Unit: Angstrom^-1.
- Parameters:
initial (bool) – Read the initial or last lattice. Useful in case of geometry optimization.
- Returns:
self.primitive_reciprocal_lattice (np.ndarray) – Lattice of the system.
- get_config_analysis(return_multiplicity=False)
Return the configuration analysis for solid solutions (CONFCON keyword in input)
- Parameters:
return_multiplicity (bool, optional) – Return multiplicity information. Defaults to False.
- Returns:
list or str – Configuration analysis if available, or a warning message
- get_convergence(history=False)
The upper level of
get_scf_convergence
andget_opt_convergence
. For analysing the geometry and energy convergence.Note
It might not work well with SCF / OPT cycles of multiple systems such as PREOPTGEOM + EOS calculations
- Parameters:
history (bool) – If true, the convergence history of optimisation (energy,gradient, displacement) / SCF is returned.
- Returns:
self (Crystal_output) – New attributes listed below
self.final_energy (float) – eV
For other attributes, see
get_scf_convergence
andget_opt_convergence
methods on the same page.
- get_scf_convergence(all_cycles=False)
Returns the scf convergence energy and energy difference. A wrapper of
CRYSTALClear.base.SCFBASE.read_convergence
.- Parameters:
all_cycles (bool, optional) – Return all SCF steps for a geometry opt. The ‘ONELOG’ CRYSTAL keyword is needed.
- Returns:
self (Crystal_output) – New attributes listed below
self.scf_cycles (int | array) – Number of cycles. Array if
all_cycles=True
.self.scf_status (str | list) – ‘terminated’, ‘converged’, ‘too many cycles’ and ‘unknown’. List if
all_cycles=True
.self.scf_energy (array) – SCF energy convergence. Unit: eV
self.scf_deltae (array) – Energy difference. Unit: eV
- get_fermi_energy(history=False)
Returns the system Fermi energy.
- Parameters:
history (bool) – Whether to read the convergence history of Fermi energy.
- Returns:
self.fermi_energy (float | array) – Fermi energy of the system. For spin-polarized insulating systems,
self.fermi_energy
would be either a 2*1 array (history=False
) or a nCYC*2 array (history=True
).
- get_band_gap(history=False)
Returns the system band gap.
- Parameters:
history (bool) – Whether to read the convergence history of band gap.
- Returns:
self.band_gap (float | array) – Band gap of the system. For spin-polarized systems,
self.band_gap
would be either a 2*1 array (history=False
) or a nCYC*2 array (history=True
).
- get_mulliken_charges()
Return the atomic Mulliken charges (PPAN keyword in input).
- Returns:
self.mulliken_charges (array) – natom*1 for non spin-polarised systems. natom*3 for spin-polarised systems. [total, \(lpha\), \(eta\)].
- get_final_energy()
Get the final energy of the system. A wrapper of
self.get_convergence
.- Returns:
self.final_energy (float) – The final energy of the system.
- get_num_cycles()
Deprecated
- Returns:
self.scf_cycles (int) – Number of SCF cycles.
- get_opt_convergence_energy()
Deprecated. Returns the energy for each opt step.
- Returns:
self.opt_energy (array) – Energy for each optimization step.
- get_opt_convergence(primitive=False, scf_history=False, write_gui=False, gui_name=None, symmetry='pymatgen', **kwargs)
Returns optimisation convergence. A wrapper of
CRYSTALClear.base.OptBASE.read_convergence
.- Parameters:
primitive (bool) – Restore the primitive cell (multiply by the cell transform matrix inversed)
scf_history (bool) – Read SCF history of each optimisation step. Keyword ‘ONELOG’ is needed. Please refer to
self.get_scf_convergence(all_cycles=True)
method.write_gui (bool) – If True, write .gui file of each step
gui_name (str) – Valid only if
write_gui = True
. Gui file is named as ‘gui_name-optxxx.gui’. If None, use ‘basename-optxxx.gui’. The basename is the same as output.symmetry (str) – Valid only if
write_gui = True
. ‘pymatgen’ to use symmetry info from a pymatgen SpacegroupAnalyzer; ‘initial’ to use symmstry information on output file. If None, no symmstry. Otherwise it is taken from the existing gui file.**kwargs – Valid only if
write_gui = True
andsymmetry = 'pymatgen'
. Passed to Pymatgen SpacegroupAnalyzer object.
- Returns:
self (Crystal_output) – New attributes listed below
self.opt_cycles (int) – Number of cycles.
self.opt_status (str) – ‘terminated’, ‘converged’, ‘failed’ and ‘unknown’
self.opt_energy (array) – Total energy convergence. Unit: eV
self.opt_deltae (array) – Total energy difference. Unit: eV
self.opt_geometry (list) – Pymatgen structure at each step.
self.opt_maxgrad (array) – Maximum gradient convergence. Unit: Hartree/Bohr
self.opt_rmsgrad (array) – RMS gradient convergence. Unit: Hartree/Bohr
self.opt_maxdisp (array) – Maximum displacement convergence. Unit: Bohr
self.opt_rmsdisp (array) – RMS displacement convergence. Unit: Bohr
- get_forces(initial=True, grad=False)
Read forces.
- Parameters:
initial (bool) – Return forces from the initial calculation. If
initial=False
, return to the last forces, which is valid only for geometry optimizations and the keyword ‘ONELOG’ is needed in d12 file.grad (bool) – Return gradient convergence history. For optimizations only.
- Returns:
self (Crystal_output) – New attributes listed below
self.forces_atoms (array) – natom*3 array. Atomic forces. Unit: Hartree/Bohr
self.forces_cell (array) – 3*3 array. Cell forces, 3D only. Unit: Hartree/Bohr
self.opt_maxgrad (array) – Maximum gradient convergence. Unit: Hartree/Bohr
self.opt_rmsgrad (array) – RMS gradient convergence. Unit: Hartree/Bohr
- get_phonon(read_eigvt=False, rm_imaginary=True, rm_overlap=True, imaginary_tol=-0.0001, q_overlap_tol=0.0001, eigvt_amplitude=1.0)
Read phonon-related properties from output file.
- Parameters:
read_eigvt (bool) – Whether to read phonon eigenvectors and normalize it to 1.
rm_imaginary (bool) – Remove the modes with negative frequencies and set all the related properties to NaN.
rm_overlap (bool) – For dispersion calculations Remove repeated q points and recalculate their weights.
imaginary_tol (float) – ``rm_imaginary`` = True only The threshold of negative frequencies.
q_overlap_tol (float) – ``rm_overlap`` = True only The threshold of overlapping points, defined as the 2nd norm of the difference of fractional q vectors
eigvt_amplitude (float | str) – ``read_eigvt = True only`` Amplitude of normalization, Or ‘classical’, ‘classical-rev’, classical amplitude and revmove classical amplitude.
Note
In QHA calculations, the ‘q point’ dimension refer to harmonic phonons computed. In other cases it refers to points in reciprocal space.
- Returns:
self (Crystal_output) – New attributes listed below
self.edft (array[float]) – \(E_{0}\) Energy with empirical correction. Unit: kJ/mol.
self.nqpoint (int) – Number of q points
self.qpoint (list[list[array[float], float]]) – A nqpoint*1 list of 2*1 list whose first element is a 3*1 array of fractional coordinates and the second is its weight.
self.nmode (array[int]) – Number of modes at q point. nqpoint*1 array.
self.frequency (array[float]) – nqpoint*nmode array ofvibrational frequency. Unit: THz
self.intens (array[float]) – nqpoint*nmode array of harmonic intensiy. Unit: km/mol
self.IR (array[bool]) – nqpoint*nmode array of boolean values specifying whether the mode is IR active
self.Raman (array[bool]) – nqpoint*nmode array of boolean values specifying whether the mode is Raman active
self.eigenvector (array[complex]) – ``read_eigvt = True only`` nqpoint*nmode*natom*3 array of eigenvectors. Normalized to 1.
- get_q_info()
Deprecated.
- get_mode()
Deprecated.
- get_phonon_eigenvector()
Deprecated.
- get_anh_spectra()
This method reads data from a CRYSTAL output file and processes it to extract anharmonic (VSCF and VCI) IR and Raman spectra.
- Returns:
self.IR_HO_0K (array[float]) – 2D numpy array containing harmonic IR frequency and intensities computed at 0 K.
self.IR_HO_T (array[float]) – 2D numpy array containing harmonic IR frequency and intensities computed at temperature T.
self.IR_VSCF_0K (array[float]) – 2D numpy array containing VSCF IR frequency and intensities computed at 0 K.
self.IR_VSCF_T (array[float]) – 2D numpy array containing VSCF IR frequency and intensities computed at temperature T.
self.IR_VCI_0K (array[float]) – 2D numpy array containing VCI IR frequency and intensities computed at 0 K.
self.IR_VCI_T (array[float]) – 2D numpy array containing VCI IR frequency and intensities computed at temperature T.
self.Ram_HO_0K_tot (array[float]) – 2D numpy array containing harmonic Raman frequency and intensities (total) computed at 0 K.
self.Ram_HO_0K_per (array[float]) – 2D numpy array containing harmonic Raman frequency and intensities (perpendicular component ) computed at temperature 0 K.
self.Ram_HO_0K_par (array[float]) – 2D numpy array containing harmonic Raman frequency and intensities (parallel component ) computed at temperature 0 K.
self.Ram_HO_T_tot (array[float]) – 2D numpy array containing harmonic Raman frequency and intensities (total) computed at temperature T.
self.Ram_HO_T_per (array[float]) – 2D numpy array containing harmonic Raman frequency and intensities (perpendicular component ) computed at temperature T.
self.Ram_HO_T_par (array[float]) – 2D numpy array containing harmonic Raman frequency and intensities (parallel component ) computed at temperature T.
self.Ram_VSCF_0K_tot (array[float]) – 2D numpy array containing VSCF Raman frequency and intensities (total) computed at 0 K. self.Ram_VSCF_0K_per (array[float]): 2D numpy array containing VSCF Raman frequency and intensities (perpendicular component) computed at 0 K.
self.Ram_VSCF_0K_par (array[float]) – 2D numpy array containing VSCF Raman frequency and intensities (parallel component) computed at 0 K.
self.Ram_VSCF_T_tot (array[float]) – 2D numpy array containing VSCF Raman frequency and intensities (total) computed at temperature T. self.Ram_VSCF_T_per (array[float]): 2D numpy array containing VSCF Raman frequency and intensities (perpendicular component) computed at temperature T.
self.Ram_VSCF_T_par (array[float]) – 2D numpy array containing VSCF Raman frequency and intensities (parallel component) computed at temperature T.
self.Ram_VCI_0K_tot (array[float]) – 2D numpy array containing VCI Raman frequency and intensities (total) computed at 0 K.
self.Ram_VCI_0K_per (array[float]) – 2D numpy array containing VCI Raman frequency and intensities (perpendicular component) computed at 0 K.
self.Ram_VCI_0K_par (array[float]) – 2D numpy array containing VCI Raman frequency and intensities (parallel component) computed at 0 K.
self.Ram_VCI_T_tot (array[float]) – 2D numpy array containing VCI Raman frequency and intensities (total) computed at temperature T.
self.Ram_VCI_T_per (array[float]) – 2D numpy array containing VCI Raman frequency and intensities (perpendicular component) computed at temperature T.
self.Ram_VCI_T_par (array[float]) – 2D numpy array containing VCI Raman frequency and intensities (parallel component) computed at temperature T.
self.Ram_HO_0K_comp_xx (array[float]) – 2D numpy array containing harmonic Raman frequency and intensities (xx component) computed at 0 K.
self.Ram_HO_T_comp_xx (array[float]) – 2D numpy array containing harmonic Raman frequency and intensities (xx component) computed at temperature T.
self.Ram_VSCF_0K_comp_xx (array[float]) – 2D numpy array containing VSCF Raman frequency and intensities (xx component) computed at 0 K.
self.Ram_VSCF_T_comp_xx (array[float]) – 2D numpy array containing VSCF Raman frequency and intensities (xx component) computed at temperature T.
self.Ram_VCI_0K_comp_xx (array[float]) – 2D numpy array containing VCI Raman frequency and intensities (xx component) computed at 0 K.
self.Ram_VCI_T_comp_xx (array[float]) – 2D numpy array containing VCI Raman frequency and intensities (xx component) computed at temperature T.
Note
Please, note that for the sake of brevity, only the xx Raman component attributes have been listed here, but the yy, zz, xy, xz, yz components are available as well.
- get_anscan(anscanwf)
Work in progress for ANSCAN stuff
- get_elatensor()
Extracts the elastic tensor from the data.
- Returns:
list – Symmetrized elastic tensor as a 6x6 nested list.
- class Properties_input(input_name=None)
Bases:
object
Create a properties_input object
- from_file(input_name)
Read the properties input from a file.
- Parameters:
input_name (str) – The name of the input file.
- Returns:
self (Properties_input) – The Properties_input object.
- make_newk_block(shrink1, shrink2, Fermi=1, print_option=0)
Returns the newk block.
- Parameters:
shrink1 (int) – The first newk shrinking factor.
shrink2 (int) – The second newk shrinking factor.
Fermi (int) – Fermi recalculation option (default is 1).
print_option (int) – Properties printing option (default is 0).
- make_bands_block(k_path, n_kpoints, first_band, last_band, print_eig=0, print_option=1, precision=5, title='BAND STRUCTURE CALCULATION')
Returns the bands block for a bands calculation.
- Parameters:
k_path (list or HighSymmKpath) – The k-path for the bands calculation.
n_kpoints (int) – The number of k-points along the path.
first_band (int) – The index of the first band.
last_band (int) – The index of the last band.
print_eig (int) – Printing options for eigenvalues (default is 0).
print_option (int) – Properties printing options (default is 1).
precision (int) – Number of zeros in the calculation of the gcd
title (str) – The title of the calculation (default is ‘BAND STRUCTURE CALCULATION’).
- make_doss_block(n_points=200, band_range=None, e_range=None, plotting_option=2, poly=12, print_option=1)
Returns the doss block for a doss calculation.
- Parameters:
n_points (int) – The number of points in the DOS plot (default is 200).
band_range (list or tuple) – The range of bands to include in the DOS calculation (default is None).
e_range (list or tuple) – The energy range for the DOS calculation (default is None).
plotting_option (int) – DOS plotting options (default is 2).
poly (int) – Degree of the polynomial for smearing (default is 12).
print_option (int) – Properties printing options (default is 1).
- make_pdoss_block(projections, proj_type='atom', output_file=None, n_points=200, band_range=None, e_range=None, plotting_option=2, poly=12, print_option=1)
Returns the pdoss block for a pdoss calculation.
- Parameters:
projections (dict) – Dictionary specifying the projections for the pdoss calculation.
proj_type (str) – Type of projection (‘atom’ or ‘site’) (default is ‘atom’).
output_file (str) – Output file name (default is None).
n_points (int) – The number of points in the DOS plot (default is 200).
band_range (list or tuple) – The range of bands to include in the DOS calculation (default is None).
e_range (list or tuple) – The energy range for the DOS calculation (default is None).
plotting_option (int) – DOS plotting options (default is 2).
poly (int) – Degree of the polynomial for smearing (default is 12).
print_option (int) – Properties printing options (default is 1).
- write_properties_input(input_name)
Writes the properties input to a file.
- Parameters:
input_name (str) – The name of the output file.
- class Properties_output
Bases:
object
Creates a Properties_output object.
- read_file(properties_output)
Parse the properties output file.
- Parameters:
properties_output (str) – The properties output file.
- Returns:
Properties_output – The updated Properties_output object.
- read_vecfield(properties_output, which_prop)
Reads the fort.25 file to return data arrays containing one or more vectiorial density properties.
- Parameters:
properties_output (str) – The properties output file.
which_prop (str) – The density property selected by the user.
'm' (magnetization), 'j' (spin current), 'J' (spin current density)
- Returns:
Properties_output (str) – The fort.25 output file.
- read_electron_band(band_file, output=None)
Generate bands object from CRYSTAL BAND.DAT or fort.25 file. Energy unit: eV. E Fermi is aligned to 0.
- Parameters:
band_file (str) – Name of BAND.DAT or fort.25 file
output (str) – Properties output file (.outp or .out). For 3D k coordinates and geometry information.
- Returns:
self.bands (BandsBASE) – A Bands base object
- read_electron_dos(properties_output)
Generate doss object from CRYSTAL DOSS.DAT or fort.25 file. Energy unit: eV. E Fermi is aligned to 0.
- Parameters:
properties_output (str) – File name
- Returns:
self.doss (DOSBASE) – A DOS base object
- read_cry_bands(properties_output)
Deprecated.
- read_cry_doss(properties_output)
Deprecated.
- read_cry_contour(properties_output)
Read the CRYSTAL contour files to create the contour objects.
- Parameters:
properties_output (str) – The properties output file.
- Returns:
Properties_output – The updated Properties_output object.
- read_cry_xrd_spec(properties_output)
Read XRD spectrum data from a file.
- Parameters:
properties_output (str) – Path to the properties output file.
- Returns:
self – The modified object with extracted XRD spectrum data.
- read_cry_rholine(properties_output)
Read density line data from a file.
- Parameters:
properties_output (str) – Path to the properties output file.
- Returns:
self – The modified object with extracted density line data.
- read_cry_seebeck(properties_output)
Read Seebeck coefficient data from a file.
- Parameters:
properties_output (str) – Path to the properties output file.
- Returns:
self – The modified object with extracted Seebeck coefficient data.
- read_cry_sigma(properties_output)
Read electrical conductivity data from a file.
- Parameters:
properties_output (str) – Path to the properties output file.
- Returns:
self – The modified object with extracted electrical conductivity data.
- read_cry_lapl_profile(properties_output)
Read Laplacian profile data from a file.
- Parameters:
properties_output (str) – Path to the properties output file.
- Returns:
self – The modified object with extracted Laplacian profile data.
- read_cry_density_profile(properties_output)
Read density profile data from a file.
- Parameters:
properties_output (str) – Path to the properties output file.
- Returns:
self – The modified object with extracted density profile data.
- read_cry_ECHG(properties_output)
Read density profile data from a file.
- Parameters:
properties_output (str) – Path to the fort.25 file.
- Returns:
self – The modified object with extracted ECHG data.
- read_cry_ECHG_delta(properties_output1, properties_output2)
Read density profile data from two files and plots the difference. It is important that the header stays the same for the two output files.
- Parameters:
properties_output1 (str) – Path to first fort.25 file.
properties_output2 (str) – Path to second fort.25 file.
- Returns:
self – The modified object with extracted ECHG data.
- class Crystal_gui
Bases:
object
This class can read a CRYSTAL gui file into an object or substrate information of the object to generate a gui file.
- read_gui(gui_file)
This is used mainly to convert the object into an ASE or pymatgen object.
- Parameters:
gui_file (str) – The CRYSTAL structure (gui) file
- write_gui(gui_file, symm=True, pseudo_atoms=[])
Write a CRYSTAL gui file (to file)
- Parameters:
gui_file (str) – The name of the gui that is going to be written ( including .gui).
symm (bool) – Whether to include symmetry operations.
pseudo_atoms (list[int]) – A list of atoms whose core is described by a pseudopotential
- class Crystal_density
Bases:
object
Read density data from a .f98 file.
- read_cry_irr_density(fort98_unit)
Read density profile data from a CRYSTAL .f98 file. :param fort98_unit: The file containing the formatted density matrix. :type fort98_unit: str
- Returns:
None
Note: This is a work in progress. If you are interested in this functionality, please open an Issue on GitHub.
- cry_combine_density(density1, density2, density3, new_density='new_density.f98', spin_pol=False)
Combine density matrix files.
- Parameters:
density1 (str) – The first density matrix file.
density2 (str) – The second density matrix file.
density3 (str) – The density matrix file for the whole system.
new_density (str, optional) – The name of the new density matrix. Defaults to ‘new_density.f98’.
spin_pol (bool, optional) – Specifies if the system is spin polarized. Defaults to False.
- Returns:
None
Note
This is a work in progress. If you are interested in this functionality, please open an Issue on GitHub.
- write_cry_density(fort98_name, new_p, new_fort98)
Write the formatted density matrix.
- Parameters:
fort98_name (str) – The name of the previous density matrix file.
new_p (list) – The new density matrix.
new_fort98 (str) – The name of the new density matrix file.
Returns
None
Note
This is a work in progress. If you are interested in this functionality, please open an Issue on GitHub.
- class External_unit
Bases:
object
- read_external_unit(external_unit)
- read_cry_irspec(external_unit)
Reads the IRSPEC.DAT unit produced by CRYSTAL ir spectra computation
- Parameters:
external_unit (str) – path to the IRSPEC.DAT file
- Returns:
The crystal_object necessary to plot the simulated spectra with plot_cry_irspec()
- read_cry_ramspec(external_unit)
Reads the RAMSPEC.DAT unit produced by CRYSTAL raman spectra computation
- Parameters:
external_unit (str) – path to the RAMSPEC.DAT file
- Returns:
The crystal_object necessary to plot the simulated spectra with plot_cry_ramspec()
- read_phonon_band(band_file, output=None)
Generate bands object from CRYSTAL PHONBAND.DAT or fort.25 file. Energy unit: eV. E Fermi is aligned to 0.
- Parameters:
band_file (str) – Name of PHONBAND.DAT or fort.25 file
output (str) – Properties output file (.outp or .out). For 3D k coordinates and geometry information.
- Returns:
self.bands (BandsBASE) – A Bands base object
- read_phonon_dos(external_unit)
Generate doss object from CRYSTAL PHONONDOSS.DAT or fort.25 file. Energy unit: eV. E Fermi is aligned to 0.
- Parameters:
properties_output (str) – File name
- Returns:
self.doss (DOSBASE) – A DOS base object
CRYSTALClear.execute
- set_runcry_path(path)
Set the path for the Runcry executable.
- Parameters:
path (str) – The path to the Runcry executable.
- Returns:
None
- set_runprop_path(path)
Set the path for the Runprop executable.
- Parameters:
path (str) – The path to the Runprop executable.
- Returns:
None
- runcry(file_name, guessp=None)
Run Runcry calculation.
- Parameters:
file_name (str) – The name of the file to run the calculation.
guessp (str, optional) – The guessp parameter. Default is None.
- Returns:
str – The result of the calculation or an error message.
- runprop(prop_name, wf_file)
Run Runprop calculation.
- Parameters:
prop_name (str) – The name of the property to calculate.
wf_file (str) – The name of the wavefunction file.
- Returns:
str – The result of the calculation or an error message.
CRYSTALClear.flatmaps
CRYSTALClear.geometry
Basic methods to manipulate pymatgen geometries
- rotate_lattice(struc, rot)
Geometries generated by
base.crysout.GeomBASE
might be rotated. Rotate them back to make them consistent with geometries in output.\[\mathbf{L}_{crys} = \mathbf{L}_{pmg}\mathbf{R}\]\(\mathbf{R}\) is the rotation matrix.
- Parameters:
struc (Structure) – Pymatgen structure
rot (array) – 3*3 numpy array, rotation matrix.
- Returns:
struc_new (Structure) – Pymatgen structure
- refine_geometry(struc, **kwargs)
Get refined geometry. Useful when reducing the cell to the irrducible one. 3D only.
- Parameters:
struc (Structure) – Pymatgen structure
**kwargs –
Passed to Pymatgen SpacegroupAnalyzer object.
- Returns:
sg (int) – Space group number
struc5 (Structure) – Irrducible structure that is consistent with International Crystallographic Table
latt (list) – minimal set of crystallographic cell parameters
natom_irr (int) – number of irrducible atoms
atom (list) – natom*4 array. 1st element: atomic number; 2-4: fractional coordinates
- get_pcel(struc, smx)
Restore the supercell to primitive cell, with the origin shifted to the middle of lattice to utilize symmetry (as the default of CRYSTAL).
- Parameters:
struc (Structure) – Pymatgen structure of supercell
smx (array) – 3*3 array of supercell expansion matrix
- Returns:
pcel (Structure) – Pymatgen structure of primitive cell
- get_scel(struc, smx)
Get the supercell from primitive cell, with the origin shifted to the middle of lattice to utilize symmetry (as the default of CRYSTAL).
- Parameters:
struc (Structure) – Pymatgen structure of supercell
smx (array) – 3*3 array of supercell expansion matrix
- Returns:
scel (Structure) – Pymatgen structure of supercell
- get_sg_symmops(struc, **kwargs)
Get space group number and corresponding symmetry operations. To keep consistency with International Crystallographic Table, refined geometry is suggested.
- Parameters:
struc (Structure) – Pymatgen Structure object.
**kwargs – Passed to Pymatgen SpacegroupAnalyzer object.
- Returns:
sg (int) – Space group number
n_symmops (int) – number of symmetry operations
symmops (array) – n_symmops*4*3 array of symmetry operations
- Miller_norm(struc, miller, d=1.0)
Find the norm vector of a specified Miller plane
- Parameters:
struc (Structure) – Pymatgen Structure object.
miller (array | list) – 3*1 list of Miller index
d (fload) – Length of norm vector
- Returns:
vec (array) – 3*1 norm vector, normalized to 1.
CRYSTALClear.plot
Functions to visualize physical properties computed with CRYSTAL .
- plot_dens_ECHG(obj_echg, levels=150, xticks=5, yticks=5, cmap_max=None, cmap_min=None, dpi=400, savefig=False, name='echg_map')
Plots the 2D ECHG density map from a fort.25 file.
- Args:
obj_echg (crystal_io.Properties_output): Properties output object. levels (int or array-like, optional): Determines the number and positions of the contour lines/regions. Default is 150. xticks (int, optional): Number of ticks in the x direction. Default is 5. yticks (int, optional): Number of ticks in the y direction. Default is 5. cmap_max(float, optional): Maximum value used for the colormap. Default is None. cmap_min(float, optional): Minimun value used for the colormap. Default is None. dpi (int, optional): DPI (dots per inch) for the output image. Default is 400. savefig (bool): Chose to save the figure or not. Default is False. name (str): Name for the colormap.
- Returns:
None
- plot_vecfield2D_m(header, dens, quivscale, name='MAG', levels=150, dpi=400)
Plots the 2D magnetization vector field.
- Parameters:
header (list) – List containing information about the fort.25 header.
dens (numpy.ndarray) – Array containing the vector field data.
quivscale (float) – Scale factor for the quiver plot.
name (str, optional) – Name used for saving the plots.
levels (int or array-like, optional) – Determines the number and positions of the contour lines/regions.
dpi (int, optional) – DPI (dots per inch) for the output image. Default is 400.
- Returns:
None
- plot_vecfield2D_j(header, dens, quivscale, name='SC', levels=150, dpi=400)
Plots the 2D vector field of the spin current.
- Parameters:
header (list) – List containing information about the fort.25 header.
dens (numpy.ndarray) – Array representing the vector field.
quivscale (float) – Scale factor for the quiver plot.
name (str, optional) – Name used for saving the plots.
levels (int or array-like, optional) – Determines the number and positions of the contour lines/regions.
dpi (int, optional) – DPI (dots per inch) for the output image. Defaults to 400.
- Returns:
None
- plot_vecfield2D_J(header, dens_JX, dens_JY, dens_JZ, quivscale, name='SCD', levels=150, dpi=400)
Plots the 2D spin current density vector fields.
- Parameters:
header (list) – List containing information about the fort.25 header.
dens_JX (numpy.ndarray) – Array representing the X-component of the spin current density.
dens_JY (numpy.ndarray) – Array representing the Y-component of the spin current density.
dens_JZ (numpy.ndarray) – Array representing the Z-component of the spin current density.
quivscale – Scale factor for the quiver plot.
name (str, optional) – Name used for saving the plots.
levels (int or array-like, optional) – Determines the number and positions of the contour lines/regions.
dpi (int, optional) – DPI (Dots per inch) for saving the plots. Defaults to 400.
- Returns:
None
- plot_phonon_band(bands, unit='cm-1', k_labels=None, mode='single', not_scaled=False, freq_range=None, k_range=None, color='blue', labels=None, linestl='-', linewidth=1, line_freq0=None, title=None, figsize=None, scheme=None, sharex=True, sharey=True, fontsize=12)
A wrapper of plot_cry_bands for phonon band structure.
- Parameters:
bands (BandsBASE|list) – Bands object generated by CRYSTALClear.crystal_io.Crystal_output.read_pband or a list of BandsBASE objects.
unit (str) – The unit of frequency. Can be ‘cm-1’ or ‘THz’.
k_labels (list) – A list of high-symmetric k point labels. Greek alphabets should be, for example, ‘Gamma’.
mode (str) – The plotting mode. Possible values are ‘single’, ‘multi’, and ‘compare’.
not_scaled (bool) – Whether to scale the x-axis for different volumes.
energy_range (array) – A 2x1 array specifying the energy range.
k_range (array) – A 2x1 array specifying the k-range.
color (str|list) – Color of plot lines. Should be consistent with bands.
labels (str|list) – Plot legend. Should be consistent with bands.
linestl (str|list) – Linestyle string. Should be consistent with bands.
linewidth (float) – The width of the plot lines.
line_freq0 (str) – The color of the frequency=0 line.
title (str) – The title of the plot.
figsize (list) – The figure size specified as [width, height].
scheme (list|tuple) – The layout of subplots.
sharex (bool) – Whether to share the x-axis among subplots.
sharey (bool) – Whether to share the y-axis among subplots.
dpi (int) – Dots per inch resolution of the saved file.
fontsize (int) – Fontsize of the axis labels.
transparency (bool) – Background transparency of the saved file,
- Returns:
Matplotlib object
- Raises:
ValueError – If the specified unit is unknown.
- plot_electron_band(bands, unit='eV', k_labels=None, mode='single', not_scaled=False, energy_range=None, k_range=None, color='blue', labels=None, linestl='-', linewidth=1, fermi='forestgreen', fermiwidth=1.5, fermialpha=1, title=None, figsize=None, scheme=None, sharex=True, sharey=True, fontsize=12)
A wrapper of plot_cry_bands for electron band structure.
- Parameters:
bands (BandsBASE|list) – Bands object generated by CRYSTALClear.crystal_io.Properties_output.read_bands or a list of BandsBASE objects.
unit (str) – The unit of energy. Can be ‘eV’ or ‘Hartree’.
k_labels (list) – A list of high-symmetric k point labels. Greek alphabets should be, for example, ‘Gamma’.
mode (str) – The plotting mode. Possible values are ‘single’, ‘multi’, and ‘compare’.
not_scaled (bool) – Whether to scale the x-axis for different volumes.
energy_range (array) – A 2x1 array specifying the energy range.
k_range (array) – A 2x1 array specifying the k-range.
color (str|list) – Color of plot lines. Should be consistent with bands.
labels (str|list) – Plot legend. Should be consistent with bands.
linestl (str|list) – Linestyle string. Should be consistent with bands.
linewidth (float) – The width of the plot lines.
fermi (str) – The color of the Fermi level line.
fermiwidth (float) – The width of the fermi line.
fermialpha (float) – Opacity of the fermi level 0-1.
title (str) – The title of the plot.
figsize (list) – The figure size specified as [width, height].
scheme (list|tuple) – The layout of subplots.
sharex (bool) – Whether to share the x-axis among subplots.
sharey (bool) – Whether to share the y-axis among subplots.
dpi (int) – Dots per inch resolution of the saved file.
fontsize (int) – Fontsize of the axis labels
transparency – Background Transparency of the saved file.
- Returns:
Matplolib object
- Raises:
ValueError – If the specified unit is unknown.
- plot_electron_dos(doss, unit='eV', beta='up', overlap=False, prj=None, energy_range=None, dos_range=None, color='blue', labels=None, linestl=None, linewidth=1, fermi='forestgreen', title=None, figsize=None)
A wrapper of plot_cry_doss for electron density of states.
- Parameters:
doss (DOSBASE) – DOS obect generated by code:CRYSTALClear.crystal_io.Properties_output.read_doss. Or a list of DOSBASE objects.
unit (str) – ‘eV’ or ‘Hartree’
beta (str) – Plot spin-down state ‘up’ or ‘down’
overlap (bool) – Plotting multiple lines into the same figure
prj (list) – Index of selected projection. Consistent with the index of the 2nd dimension of
doss.doss
energy_range (list[float]) – 2*1 list of energy range
dos_range (list[float]) – 2*1 list of DOS range
color (str | list[str]) – Color of plot lines. Should be consistent with number of projections.
labels (str | list[str]) – Plot legend. Should be consistent with number of projections.
linestl (str | list[str]) – linestyle string. Should be consistent with number of projections.
linewidth (float)
fermi (str) – Color of Fermi level line.
title (str)
figsize (list[float])
- Returns:
Matplotlib object
- plot_phonon_dos(doss, unit='cm-1', overlap=False, prj=None, freq_range=None, dos_range=None, color='blue', labels=None, linestl=None, linewidth=1, line_freq0=None, title=None, figsize=None)
A wrapper of plot_cry_doss for electron density of states.
- Parameters:
doss (DOSBASE) – DOS obect generated by code:CRYSTALClear.crystal_io.Crystal_output.read_pdos. Or a list of DOSBASE objects.
unit (str) – ‘cm-1’ or ‘THz’
overlap (bool) – Plotting multiple lines into the same figure
prj (list) – Index of selected projection. Consistent with the index of the 2nd dimension of
doss.doss
freq_range (list[float]) – 2*1 list of frequency range
dos_range (list[float]) – 2*1 list of DOS range
color (str | list[str]) – Color of plot lines. Should be consistent with number of projections.
labels (str | list[str]) – Plot legend. Should be consistent with number of projections.
linestl (str | list[str]) – linestyle string. Should be consistent with number of projections.
linewidth (float)
line_freq0 (str) – Color of frequency = 0 line.
title (str)
figsize (list[float])
- Returns:
Matplotlib object
- plot_electron_banddos(bands, doss, unit='eV', k_labels=None, dos_beta='down', dos_prj=None, energy_range=None, dos_range=None, color_band='blue', color_dos='blue', labels=None, linestl_band='-', linestl_dos=None, linewidth=1, fermi='forestgreen', title=None, figsize=None)
A wrapper of plot_cry_es for electron band structure + dos. For spin-polarized cases, beta state.
- Parameters:
bands (BandsBASE|list) – Bands object generated by CRYSTALClear.crystal_io.Properties_output.read_bands or a list of BandsBASE objects.
doss (DOSBASE) – DOS object generated by CRYSTALClear.crystal_io.Properties_output.read_doss or a list of DOSBASE objects.
unit (str) – Unit of energy. Valid options are ‘eV’ or ‘Hartree’.
k_labels (list) – A list of high-symmetric k point labels. Greek alphabets should be represented as strings, for example, ‘Gamma’.
dos_beta (str) – Spin state to plot. Valid options are ‘Up’ or ‘down’. If ‘down’, the beta state will be plotted on the same side as the alpha state, otherwise on the other side.
dos_prj (list) – Index of selected projection. Consistent with the index of the 2nd dimension of doss.doss.
energy_range (list) – A list of two values representing the energy range to be plotted.
dos_range (list) – DOS range for the y-axis.
color_band (str) – Color of the electron bands in the plot.
color_dos (str) – Color of the density of states (DOS) in the plot.
labels (list) – A list of labels for the plot legend.
linestl_band (str) – Linestyle of the electron bands.
linestl_dos (str) – Linestyle of the density of states (DOS).
linewidth (float) – Width of the lines in the plot.
fermi (str) – Color of the Fermi level line.
title (str) – Title of the plot.
figsize (list[float]) – Size of the figure in inches (width, height).
legend (bool) – Enables or disables the legend of the density of states (DOS).
- Returns:
Matplotlib object
- Raises:
ValueError – If the unit parameter is unknown.
- plot_phonon_banddos(bands, doss, unit='cm-1', k_labels=None, dos_prj=None, freq_range=None, dos_max_range=None, color_band='blue', color_dos='blue', labels=None, linestl_band='-', linestl_dos=None, linewidth=1, freq0_line=None, title=None, figsize=None)
A wrapper of plot_cry_es for phonon band structure + dos. Only one pair is permitted.
- Parameters:
bands (BandsBASE|list) – Bands object generated by CRYSTALClear.crystal_io.Properties_output.read_bands or a list of BandsBASE objects.
doss (DOSBASE) – DOS object generated by CRYSTALClear.crystal_io.Properties_output.read_doss or a list of DOSBASE objects.
unit (str) – Unit of frequency. Valid options are ‘cm-1’ or ‘THz’.
k_labels (list) – A list of high-symmetric k point labels. Greek alphabets should be represented as strings, for example, ‘Gamma’.
dos_prj (list) – Index of selected projection. Consistent with the index of the 2nd dimension of doss.doss.
freq_range (list) – A list of two values representing the frequency range to be plotted.
dos_max_range (float) – Maximum DOS range for the y-axis.
color_band (str) – Color of the phonon bands in the plot.
color_dos (str) – Color of the density of states (DOS) in the plot.
labels (list) – A list of labels for the plot legend.
linestl_band (str) – Linestyle of the phonon bands.
linestl_dos (str) – Linestyle of the density of states (DOS).
linewidth (float) – Width of the lines in the plot.
freq0_line (str) – Color of the frequency=0 line.
title (str) – Title of the plot.
figsize (list[float]) – Size of the figure in inches (width, height).
- Returns:
Matplotlib object
- Raises:
ValueError – If the unit parameter is unknown.
- plot_cry_contour(contour_obj)
Plot a contour plot.
- Parameters:
contour_obj (object) – Contour object representing the contour plot.
- Returns:
None
Notes
Plots a contour plot based on the data in the contour object.
Retrieves the data from the contour object and converts it to a 2D list.
Sets the figure size based on x_graph_param and y_graph_param attributes of the contour object.
Sets the x-axis and y-axis labels.
Creates a meshgrid using the x_points and y_points attributes of the contour object.
Defines contour levels, colors, linestyles, and fmt.
Plots the contour plot.
Saves the plot to a file named ‘figure_TIPO_YYYY-MM-DD_HHMMSS.jpg’ in the current directory.
- plot_cry_contour_differences(contour_obj, contour_obj_ref)
Plot the differences between two contour plots.
- Parameters:
contour_obj (object) – Contour object representing the original contour plot.
contour_obj_ref (object) – Contour object representing the reference contour plot.
- Returns:
None
Notes
Plots the differences between two contour plots.
Requires the contour objects to have a tipo attribute with values ‘SURFLAPP’, ‘SURFLAPM’, ‘SURFRHOO’, or ‘SURFELFB’.
Calculates the difference between the dataframes of the two contour objects.
Sets the figure size based on x_graph_param and y_graph_param attributes of the contour object.
Sets the x-axis and y-axis labels.
Creates a meshgrid using the x_points and y_points attributes of the contour object.
Defines contour levels, colors, and linestyles.
Plots the contour differences.
Saves the plot to a file named ‘figure_diff_TIPO_YYYY-MM-DD_HHMMSS.jpg’ in the current directory.
- plot_cry_xrd(xrd_obj)
Plot the X-ray diffraction pattern.
- Parameters:
xrd_obj (object) – XRD object containing the data for the X-ray diffraction pattern.
- Returns:
None
Notes
Plots the X-ray diffraction pattern.
Sets the figure size to [16, 9].
Sets the x-axis limit to (0, 30).
Saves the plot to a file named ‘figure_XRD_YYYY-MM-DD_HHMMSS.jpg’ in the current directory.
- plot_cry_rholine(rholine_obj)
Plot the resistivity as a function of distance.
- Args:
rholine_obj (object): Rholine object containing the data for the resistivity.
- Returns:
None
- Notes:
Plots the resistivity as a function of distance.
Sets the x-axis label as ‘d [$AA$]’ and the y-axis label as r’$
ho$ [$ rac{e}{AA^3}$]’.
Saves the plot to a file named ‘figure_rholine_YYYY-MM-DD_HHMMSS.jpg’ in the current directory.
- plot_cry_lapl_profile(lapl_obj)
Plot the Laplacian profile of a crystal.
- Parameters:
lapl_obj (object) – Laplacian object containing the data for the Laplacian profile.
- Returns:
None
Notes
Plots the Laplacian profile using the data from the Laplacian object.
The x-axis represents the distance in angstroms.
The y-axis represents the Laplacian in electrons per cubic angstrom to the fifth power (e/A^5).
The area under the curve where the Laplacian is negative is filled with a light blue color.
The area under the curve where the Laplacian is positive is filled with a light coral color.
- plot_cry_density_profile(lapl_obj)
Plot the density profile of a crystal.
- Parameters:
lapl_obj (object) – Laplacian object containing the data for the density profile.
- Returns:
None
Notes
Plots the density profile using the data from the Laplacian object.
The x-axis represents the distance in angstroms.
The y-axis represents the density in electrons per cubic angstrom (e/A^3).
- plot_cry_seebeck_potential(seebeck_obj)
Plot the Seebeck coefficient as a function of chemical potential.
- Parameters:
seebeck_obj (object) – Seebeck object containing the data for the Seebeck coefficient.
- Returns:
None
Notes
Prompts the user to choose the direction to plot among S_xx, S_xy, S_xz, S_yx, S_yy, S_yz, S_yz, S_zx, S_zy, S_zz.
Plots the Seebeck coefficient as a function of chemical potential for each temperature.
Distinguishes between n-type and p-type conduction with dashed and solid lines, respectively.
- plot_cry_seebeck_carrier(seebeck_obj)
Plot the Seebeck coefficient as a function of charge carrier concentration.
- Parameters:
seebeck_obj – Seebeck object containing the data for the Seebeck coefficient.
- Returns:
None
Notes
Prompts the user to choose the direction to plot among S_xx, S_xy, S_xz, S_yx, S_yy, S_yz, S_yz, S_zx, S_zy, S_zz.
Plots the Seebeck coefficient as a function of charge carrier concentration for each temperature, distinguishing between n-type and p-type conduction.
- plot_cry_multiseebeck(*seebeck)
Plot the multiseebeck coefficient for different temperatures.
- Parameters:
*seebeck – Variable number of seebeck objects containing the data for the Seebeck coefficient.
- Returns:
None
Notes
Prompts the user to input the index of the temperature to plot.
Prompts the user to input the lower and higher values of chemical potential to plot in eV.
Prompts the user to choose the direction to plot among S_xx, S_xy, S_xz, S_yx, S_yy, S_yz, S_yz, S_zx, S_zy, S_zz.
Plots the multiseebeck coefficient for each seebeck object.
Differentiates transport coefficients due to n-type or p-type conduction using dashed and solid lines.
Saves the plot to a file named ‘multiseebeckYYYY-MM-DD_HHMMSS.jpg’, where YYYY-MM-DD_HHMMSS represents the current date and time.
- plot_cry_sigma_potential(sigma_obj)
Plot the electrical conductivity as a function of chemical potential.
- Parameters:
sigma_obj (object) – Sigma object containing the data for electrical conductivity.
- Returns:
None
Notes
Prompts the user to choose the direction to plot among S_xx, S_xy, S_xz, S_yy, S_yz, S_zz.
Plots the electrical conductivity as a function of chemical potential for each temperature.
Distinguishes between n-type and p-type conduction with dashed and solid lines, respectively.
- plot_cry_sigma_carrier(sigma_obj)
Plot the electrical conductivity as a function of charge carrier concentration.
- Parameters:
sigma_obj – Sigma object containing the data for the electrical conductivity.
- Returns:
None
Notes
Prompts the user to choose the direction to plot among S_xx, S_xy, S_xz, S_yy, S_yz, S_zz.
Plots the electrical conductivity as a function of charge carrier concentration for each temperature, distinguishing between n-type and p-type conduction.
- plot_cry_multisigma(*sigma)
Plot the multisigma conductivity for different temperatures.
- Parameters:
*sigma – Variable number of sigma objects containing the data for the conductivity.
- Returns:
None
Notes
Prompts the user to input the index of the temperature to plot.
Prompts the user to input the lower and higher values of chemical potential to plot in eV.
Prompts the user to choose the direction to plot among S_xx, S_xy, S_xz, S_yy, S_yz, S_zz.
Plots the multisigma conductivity for each sigma object.
Differentiates transport coefficients due to n-type or p-type conduction using dashed and solid lines.
Saves the plot to a file named ‘multisigmaYYYY-MM-DD_HHMMSS.jpg’, where YYYY-MM-DD_HHMMSS represents the current date and time.
- plot_cry_powerfactor_potential(seebeck_obj, sigma_obj)
Plot the power factor for different potentials.
- Parameters:
seebeck_obj – Seebeck object containing the data for the Seebeck coefficient.
sigma_obj – Sigma object containing the data for the electrical conductivity.
- Returns:
None
Notes
Prompts the user to choose the direction to plot among PF_xx, PF_xy, PF_xz, PF_yx, PF_yy, PF_yz, PF_yz, PF_zx, PF_zy, PF_zz.
Calculates the power factor using the Seebeck coefficient and electrical conductivity data for each temperature.
Plots the power factor for each temperature as a function of the chemical potential, distinguishing between n-type and p-type conduction.
- plot_cry_powerfactor_carrier(seebeck_obj, sigma_obj)
Plot the power factor for different charge carrier concentrations.
- Parameters:
seebeck_obj – Seebeck object containing the data for the Seebeck coefficient.
sigma_obj – Sigma object containing the data for the electrical conductivity.
- Returns:
None
Notes
Prompts the user to choose the direction to plot among PF_xx, PF_xy, PF_xz, PF_yx, PF_yy, PF_yz, PF_yz, PF_zx, PF_zy, PF_zz.
Calculates the power factor using the Seebeck coefficient and electrical conductivity data for each temperature.
Plots the power factor for each temperature as a function of the charge carrier concentration, distinguishing between n-type and p-type conduction.
- plot_cry_zt(seebeck_obj, sigma_obj)
Plot the ZT value for different temperatures.
- Parameters:
seebeck_obj – Seebeck object containing the data for the Seebeck coefficient.
sigma_obj – Sigma object containing the data for the electrical conductivity.
- Returns:
None
Notes
Prompts the user to input the value of ktot in W-1K-1m-1.
Prompts the user to choose the direction to plot among ZT_xx, ZT_xy, ZT_xz, ZT_yx, ZT_yy, ZT_yz, ZT_yz, ZT_zx, ZT_zy, ZT_zz.
Calculates the ZT value using the Seebeck coefficient and electrical conductivity data.
Plots the ZT value for each temperature as a function of the chemical potential.
- plot_cry_young(theta, phi, S)
Compute Young’s modulus for each direction of the space (i.e., each pair of theta and phi angles).
- Parameters:
theta (float) – Theta value.
phi (float) – Phi value.
S (numpy.ndarray) – Compliance matrix.
- Returns:
float – Young’s modulus values.
Notes
This function is intended to be called by cry_ela_plot
- plot_cry_comp(theta, phi, S)
Compute linear compressibility for each direction of the space (i.e., each pair of theta and phi angles).
- Parameters:
theta (float) – Theta value.
phi (float) – Phi value.
S (numpy.ndarray) – Compliance matrix.
- Returns:
float – Linear compressibility values.
Notes
This function is intended to be called by cry_ela_plot
- plot_cry_shear(theta_1D, phi_1D, S, ndeg, shear_choice)
For each direction of the space (i.e., for each pair of theta and phi angles) the shear modulus is computed for the third angle chi and the average, maximum and minimum values are stored.
- Parameters:
theta_1D (numpy.ndarray) – One-dimensional array of theta values.
phi_1D (numpy.ndarray) – One-dimensional array of phi values.
S (numpy.ndarray) – Compliance matrix.
ndeg (int) – Number of degrees for discretization.
shear_choice (str) – Type of shear property to plot. Options: “avg”, “min”, “max”.
- Returns:
numpy.ndarray – Shear property array.
Notes
This function is intended to be called by cry_ela_plot
- plot_cry_poisson(theta_1D, phi_1D, S, ndeg, poisson_choice)
For each direction of the space (i.e., for each pair of theta and phi angles) the Poisson ratio is computed for the third angle chi and the average, maximum and minimum values are stored.
- Parameters:
theta_1D (numpy.ndarray) – One-dimensional array of theta values.
phi_1D (numpy.ndarray) – One-dimensional array of phi values.
S (numpy.ndarray) – Compliance matrix.
ndeg (int) – Number of degrees for discretization.
poisson_choice (str) – Type of Poisson’s ratio to plot. Options: “avg”, “min”, “max”.
- Returns:
numpy.ndarray – Poisson’s ratio array.
Notes
This function is intended to be called by cry_ela_plot
- plot_cry_ela(choose, ndeg, *args)
Plot crystal elastic properties on the basis of the elastic tensor. A variable number of elastic tensors can be provided in order to get multiple plots in one shot, establishing a fixed color scale among them.
- Parameters:
choose (str) – Property to plot. Options: “young”, “comp”, “shear avg”, “shear min”, “shear max”, “poisson avg”, “poisson min”, “poisson max”.
ndeg (int) – Number of degrees for discretization.
*args – Variable number of elastic tensors.
- Returns:
Tuple of lists
- fig_list – list of matplotlib.figure.Figure A list containing matplotlib Figure objects for each plot.
- ax_list – list of matplotlib.axes._axes.Axes A list containing the Axes objects associated with each plot.
- plt_list – list of matplotlib.pyplot A list of the pyplot objects for each plot, representing the actual plot.
- plot_cry_irspec(irspec, x_unit='cm-1', y_mode='LG', figsize=None, linestyle='-', linewidth=1.5, color='tab:blue', freq_range=None, int_range=None, label=None)
Generates the IR spectra for the IRSPEC.DAT file produced by an IRSPEC calculation
- Parameters:
irspec (External_unit object) – Object generated by the read_cry_irspec function necessary for the plot
x_unit (str, optional) – Unit measure of the x axes. Avalilable: ‘cm-1’ and ‘nm’. Defaults to ‘cm-1’.
y_mode (str, optional) – Peak broadening modality in absorbance and reflectance. Available: ‘LG’(Lorentzian-Gaussian broadening), ‘V’ (Voight broadening), ‘RS’ (Rayleigh spherical particles), ‘RE’ (Rayleigh with elipsoid particles), ‘REFL’ (Reflectance) Defaults to ‘LG’.
figsize (tuple, optional) – Image dimensions correspondig to matplotlib figsize. Defaults to None.
linestyle (str/list[str], optional) – linestyle corresponding to the matplotlib one it can be a list for a multiplot. Defaults to ‘-‘.
linewidth (float/list[float], optional) – linewidth corresponding to the matplotlib one it can be a list for a multiplot. Defaults to 1.5.
color (str/list[str], optional) – Color of the spectra it can accept all matplotlib colors it can be a list for multiplots. Defaults to ‘tab:blue’.
freq_range (list, optional) – Two element list [min, max], that allows to visualize the spectra in a given frequency window. Defaults to None.
int_range (list, optional) – Two element list [min, max], that allows to visualize the spectra in a given intensity window. Defaults to None.
label (list[str], optional) – List of labels for the legend of a multiplot. Defaults to None.
dpi (int, optional) – Resolution of the saved file. Defaults to 300.
transparency (bool, optional) – Enables the transparency of the saved file background. Defaults to False.
- Returns:
Matplotlib object
- Raises:
ValueError – The function raises an error when the object to be plotted does not have the required y_mode
- plot_cry_ramspec(ramspec, y_mode='total', figsize=None, linestyle='-', linewidth=1.5, color='tab:blue', freq_range=None, int_range=None, label=None)
Generates the RAMAN spectra for the RAMSPEC.DAT file produced by an RAMSPEC calculation
- Parameters:
irspec (External_unit object) – Object generated by the read_cry_ramspec function necessary for the plot
y_mode (str, optional) – Polarization of the spectra for the simulated compound Available: ‘total’, ‘parallel’, ‘perpendicular’ (for powders), ‘xx’, ‘xy’, ‘xz’, ‘yy’, ‘yz’, ‘zz’ (for single crystals) Defaults to ‘LG’.
figsize (tuple, optional) – Image dimensions correspondig to matplotlib figsize. Defaults to None.
linestyle (str/list[str], optional) – linestyle corresponding to the matplotlib one it can be a list for a multiplot. Defaults to ‘-‘.
linewidth (float/list[float], optional) – linewidth corresponding to the matplotlib one it can be a list for a multiplot. Defaults to 1.5.
color (str/list[str], optional) – Color of the spectra it can accept all matplotlib colors it can be a list for multiplots. Defaults to ‘tab:blue’.
freq_range (list, optional) – Two element list [min, max], that allows to visualize the spectra in a given frequency window. Defaults to None.
int_range (list, optional) – Two element list [min, max], that allows to visualize the spectra in a given intensity window. Defaults to None.
label (list[str], optional) – List of labels for the legend of a multiplot. Defaults to None.
dpi (int, optional) – Resolution of the saved file. Defaults to 300.
transparency (bool, optional) – Enables the transparency of the saved file background. Defaults to False.
- Returns:
Matplotlib object
- plot_cry_spec(transitions, typeS, components=False, bwidth=5, stdev=3, eta=0.5, fmin=None, fmax=None, ylim=None, savefig=False, dpi=300, filetype='png', exp_spec=None, sep=';', show=True, export_csv=False, label=None, xlabel='Wavenumber [cm$^{-1}$]', ylabel='Intensity [arb. u.]', linewidth=2.0, padd=100, fontsize=12, style=None, compstyle=None, nopadding=False, figsize=(16, 6))
This function enables the simulation of vibrational spectra based on a 2D NumPy array containing a list of transition frequencies and the corresponding intensities. The code allows users to model spectral broadening according to various profiles (Gaussian, Lorentzian, pseudo-Voigt), or zero broadening (Dirac deltas-like lines). Please, note that by turning the optional argument ‘component’ to True you can additionally plot contributions arising from each transition.
- Parameters:
transitions (float|numpy.ndarray) – Array containing transition frequencies
intensities ((axis=0) and corresponding)
typeS (str) – String specifying the spectral profile: ‘bars’,
'lorentz'
'gauss'
'pvoigt'.
components (bool, optional) – Whether to plot contributions arising from
transition (each)
bwidth (float, optional) – Half-width at half-maximum of the Lorentzian
profile (default is 0.5)
stdev (float, optional) – Standard deviation of the Gaussian profile
5). ((default is)
eta (float, optional) – Fraction of Lorentzian character in pseudo-Voigt
profile
fmin (float, optional) – Minimum frequency.
fmax (float, optional) – Maximum frequency.
ylim (float, optional) – Maximum intensity.
savefig (bool, optional) – Whether to save the figure (default is False).
dpi (float, optional) – Dots per inches (default is 300).
filetype (str, optional) – File extension (default is ‘png’).
show (bool, optional) – Whether to show the figure (default is True).
export_csv (bool, optional) – Whether to save plot in csv format (default is
False).
xlabel (str, optional) – x-axis label (default is ‘Wavenumber [cm$^{-1}$]’).
ylabel (str, optional) – y-axis label (default is ‘Intensity [arb. u.]’).
linewidth (float) – Linewidth (default is 2.0).
padd (float, optional) – left- and right- hand side padding expressed in the
x-axis (same unit of the quantity reported in)
fontsize (integer, optional) – Fontsize (default is 12).
style (str, optional) – String specifying Matplotlib style.
compstyle (str|list, optional) – List containing Matplotlib styles to plot
component. (each)
nopadding (bool, optional) – Whether to remove padding (default is False).
figsize (real|list, optional) – List of two numbers specifying the aspect
figure (ratio of the)
- Returns:
:class:`matplotlib.pyplot`
A matplotlib object representing the result of the plot
- plot_cry_spec_multi(files, typeS, components=False, bwidth=5, stdev=3, eta=0.5, fmin=None, fmax=None, ylim=None, savefig=False, dpi=300, filetype='png', label=None, xlabel='Wavenumber [cm$^{-1}$]', ylabel='Intensity [arb. u.]', linewidth=2.0, padd=100, fontsize=12, style=None, nopadding=False, figsize=(16, 6), exp_spec=None, norm_fac=1, sep=';')
This function is a wrapper for plot_spec function, enablng the simulation of many vibrational spectra coming from a list of NumPy array.
- Parameters:
transitions (float|numpy.ndarray) – Array containing transition frequencies
intensities ((axis=0) and corresponding)
typeS (str) – String specifying the spectral profile: ‘bars’,
'lorentz'
'gauss'
'pvoigt'.
components (bool, optional) – Whether to plot contributions arising from
transition (each)
bwidth (float, optional) – Half-width at half-maximum of the Lorentzian
profile (default is 0.5)
stdev (float, optional) – Standard deviation of the Gaussian profile
5). ((default is)
eta (float, optional) – Fraction of Lorentzian character in pseudo-Voigt
profile
fmin (float, optional) – Minimum frequency.
fmax (float, optional) – Maximum frequency
ylim (float, optional) – Maximum intensity.
savefig (bool, optional) – Whether to save the figure (default is False).
dpi (float, optional) – Dots per inches (default is 300).
filetype (str, optional) – File extension (default is ‘png’).
xlabel (str, optional) – x-axis label (default is ‘Wavenumber [cm$^{-1}$]’).
ylabel (str, optional) – y-axis label (default is ‘Intensity [arb. u.]’).
linewidth (float) – Linewidth (default is 2.0).
padd (float, optional) – left- and right- hand side padding expressed in the
x-axis (same unit of the quantity reported in)
fontsize (integer, optional) – Fontsize (default is 12).
style (str, optional) – String specifying Matplotlib style.
nopadding (bool, optional) – Whether to remove padding (default is False).
figsize (real|list, optional) – List of two numbers specifying the aspect
figure (ratio of the)
- Returns:
:class:`matplotlib.pyplot`
A matplotlib object representing the result of the plot
- plot_cry_anscan(co, scale_wf=None, scale_prob=None, harmpot=False, scanpot=True, figsize=[10, 10])
This function provides a plotting tool for the ANSCAN keyword.
- Parameters:
co (crystal_io.Crystal_output) – Crystal output object.
scale_wf (float, optional) – Scaling factor for wavefunctions plot. By
default (plot. By)
this (wavefunctions are not plotted. Providing a value for)
accordingly. (scales it)
scale_prob (float, optional) – Scaling factor for probability density
default
a (probability densities are not plotted. Providing)
and (value for this argument enables the probability density plot)
accordingly.
harmpot (bool, optional) – A logical flag to activate the plotting of
potential (the harmonic)
scanpot (bool, optional) – A logical flag to activate the plotting of
ANSCAN (the scan potential provided by)
figsize (list[float])
- Returns:
:class:`matplotlib.pyplot`
A matplotlib object representing the result of the plot
CRYSTALClear.thermodynamics
A post-processing module for DFT lattice dynamics by harmonic and quasiharmonic approximations (HA/QHA).
- class Mode(rank=0, frequency=[], volume=[], eigenvector=[])
Bases:
object
Defines a vibrational mode and do analysis per mode.
- Parameters:
rank (int) – The rank of the mode object, from 1.
frequency (array[float] | list[float]) – Frequencies of the mode (Ncalc*1). Unit: THz. Note: NOT angular frequency, which is frequency * 2pi.
volume (array[float] | list[float]) – Lattice volumes of harmonic calculations (Ncalc*1). Unit: Angstrom^3
eigenvector (array[float] | list[float]) – Corresponding normalized eigenvectors (Ncalc*Natom*3).
- Returns:
self.rank (int)
self.ncalc (int) – The number of harmonic calculations (typically at different volumes)
self.frequency (array[float])
self.volume (array[float])
self.eigenvector (array[float])
- get_zp_energy()
Get the zero-point energy of a single mode. ncalc = 1 only.
\[E^{zp}_{i,\mathbf{q}}=\frac{1}{2}\hbar\omega_{i,\mathbf{q}}\]- Returns:
self.zp_energy (float) – Zero-point energy. Unit: KJ/mol
- get_u_vib(temperature=298.15)
Get the vibration contribution to internal energy (including zero-point energy) of a single mode. ncalc = 1 only.
\[U^{vib}_{i,\mathbf{q}}\left(T\right)=E^{zp}_{i,\mathbf{q}}+ \frac{\hbar\omega_{i,\mathbf{q}}}{\exp{\left( \frac{\hbar\omega_{i,\mathbf{q}}}{k_{B}T} \right)}-1}\]- Parameters:
temperature (float, optional) – Temperature where the quantity is computed. Unit: K
- Returns:
self.u_vib (float) – Vibration contribution to internal energy. Unit: KJ/mol
- get_entropy(temperature)
Get the entropy of a single mode. ncalc = 1 only.
\[S_{i,\mathbf{q}}\left(T\right)=k_{B}\left\{ \frac{\hbar\omega_{i,\mathbf{q}}}{k_{B}T\left[ \exp{\left( \frac{\hbar\omega_{i,\mathbf{q}}}{k_{B}T} \right)}-1 \right]}-\ln{\left[ 1-\exp{\left( -\frac{\hbar\omega_{i,\mathbf{q}}}{k_{B}T} \right)} \right]} \right\}\]- Parameters:
temperature (float, optional) – Unit: K
- Returns:
self.entropy (float) – Entropy. Unit: J/mol*K
- get_c_v(temperature)
Get the constant volume specific heat of a single mode. ncalc = 1 only.
\[C^{V}_{i,\mathbf{q}}= \frac{\left(\hbar\omega_{i,\mathbf{q}}\right)^{2}}{k_{B}T^{2}} \frac{\exp{ \left( \frac{\hbar\omega_{i,\mathbf{q}}}{k_{B}T} \right)} }{\left[ \exp{\left( \frac{\hbar\omega_{i,\mathbf{q}}}{k_{B}T} \right)-1} \right]^{2} }\]- Parameters:
temperature (float, optional) – Unit: K
- Returns:
self.c_v (float) – Constant volume specific heat. Unit: J/mol*K
- get_classical_amplitude(struc)
Empty method - development ongoing
- polynomial_fit(order=[2, 3])
Fit phonon frequency as the polynomial function of volume. ncalc > 1 only.
Note
To improve the accuracy of fittings, \(\Delta\omega(\Delta V)\) is fitted as a polynomial function without the constant term.
\(\Delta V=V-V_{min}\) is used so HA phonons of the most compact structure is kept. See ‘FIXINDEX’ keyword in CRYSTAL manual for further information.
- Parameters:
order (array[int] | list[int]], optional) – Orders of polynomials.
- Returns:
self.poly_fit (Dict[int, NumPy Polynomial]) – Key - orders of power, Value - fitted NumPy polynomials
self.poly_fit_rsquare (Dict[int, float]) – Key - orders of power, Value - goodness of fittings, characterized by R^2.
- get_gruneisen(order, volume)
Return to mode Grüneisen parameter. ncalc > 1 only.
\[\gamma = -\frac{V}{\omega(V)}\frac{\partial\omega}{\partial V}\]- Parameters:
order (int | list[int]) – See
polynomial_fit
volume (float | array) – Typically the equilibrium volume
- Returns:
self.gruneisen (dict) – Key, order; Value, Gruneisen parameter
- class Harmonic(temperature=[], pressure=[], filename=None, autocalc=True)
Bases:
object
A class for harmonic phonon calclulations. It can be parameterized from a CRYSTAL output file, phonopy ouput file or by setting all the information (usually for QHA).
- Parameters:
temperature (array[float] | list[float], optional) – Temperatures where thermodynamic properties are computed. Unit: K
pressure (array[float] | list[float], optional) – Pressures where the thermodyanmic properties are calculated. Unit: GPa
filename (str | None) – Name of the printed-out file. If None, do not print out file.
autocalc (bool) – Automatically launch calculations.
Temperatures and pressures can also be defined by
self.thermodynamics
, whose entries always cover the entries here.Phonon dispersions are forced to be summed if the automatic scheme (
autocalc=True
) is launched. To get verbose outputs, callself.thermodynamics()
first and then callself.print_results()
.Usage:
ha = Harmonic(temperature=[0, 100, 200, 300], pressure=[0.,]) ha.from_file('harmonic_phonon.out')
- from_file(output_name, scelphono=[], read_eigvt=False, imaginary_tol=-0.0001, q_overlap_tol=0.0001)
Generate the Harominc object from a HA output file. Imaginary modes and overlapped q points are forced to be cleaned.
- Parameters:
output_name (str) – Name of the output file.
scellphono (array[float] | list[float], optional) – The ‘SCELPHONO’ keyword in CRYSTAL input file. By default a 1*1*1 ‘SCELPHONO’ is assumed.
read_eigvt (bool) – Whether to read eigenvectors from output.
imaginary_tol (float) – The threshold of negative frequencies.
q_overlap_tol (float) – The threshold of overlapping points, defined as the 2nd norm of the difference of fractional q vectors
- Returns:
self.structure (PyMatGen Structure) – Cell reduced by SCELPHONO.
self.natom (int) – Number of atoms in the reduced cell.
self.volume (float) – Volume of the reduced cell. Unit: Angstrom^3
self.edft (float)
self.nqpoint (int)
self.qpoint (list)
self.nmode (array[int])
self.mode (list[Mode]) – List of mode objects at all the qpoints.
- Raises:
ValueError – If a QHA output file is read.
- from_phonopy(phono_yaml, struc_yaml=None, edft=None, scale=1.0, imaginary_tol=-0.0001, q_overlap_tol=0.0001, q_id=None, q_coord=None)
Build a Harmonic object from Phonopy ‘band.yaml’ or ‘qpoints.yaml’ file.
- Parameters:
phono_yaml (str) – Phonopy band.yaml or qpoint.yaml file
struc_yaml (str) – Phonopy phonopy.yaml or phonopy_disp.yaml file. Needed only if a qpoint.yaml file is read.
edft (float) – DFT energy. Unit: kJ/mol
scale (float) – Scaling factor of phonon frequency.
imaginary_tol (float) – The threshold of negative frequencies.
q_overlap_tol (float) – The threshold of overlapping points, defined as the 2nd norm of the difference of fractional q vectors
q_id (list[int]) – Specify the id (from 0) of q points to be read. nqpoint*1 list.
q_coord (list[list]) – Specify the coordinates of q points to be read. nqpoint*3 list.
q_id
andq_coord
should not be set simultaneously. If set,q_id
takes priority andq_coord
is ignored. If both are none, all the points will be read.- Raises:
Exception – If the length unit in yaml file is neither ‘au’ nor ‘angstrom’.
Exception – If q point is not found.
- from_frequency(edft, qpoint, frequency, eigenvector, structure=None, natom=None, volume=None, imaginary_tol=-0.0001, q_overlap_tol=0.0001, ignore_natom=False)
Generate a Harmonic object by specifying frequency and eigenvector. Imaginary modes and overlapped q points are forced to be cleaned.
- Parameters:
edft (float) – Electron total energy
qpoint (list[list[array[float], float]]) – Fractional coordinate and weight of qpoint
frequency (array[float]) – Array of frequencies. Unit: THz
eigenvector (array[float]) – Normalized eigenvectors.
structure (Pymatgen Structure)
natom (int)
volume (float)
imaginary_tol (float) – The threshold of negative frequencies.
q_overlap_tol (float) – The threshold of overlapping points, defined as the 2nd norm of the difference of fractional q vectors
ignore_natom (bool) – Developer only.
Note
The user should define either
structure
ornatom
+volume
.- Returns:
self.structure (PyMatGen Structure) – Cell reduced by SCELPHONO.
self.natom (int) – Number of atoms in the reduced cell.
self.volume (float) – Volume of the reduced cell. Unit: Angstrom^3
self.edft (float)
self.nqpoint (int)
self.qpoint (list)
self.nmode (array[int])
self.mode (list[Mode]) – List of mode objects at all the qpoints.
- Raises:
AttributeError – If computational data is stored in the object.
ValueError – If neither of the 2 available options are defined.
- thermodynamics(sumphonon=True, mutewarning=False, **kwargs)
Calculate the thermodynamic properties (zp_energy, u_vib, entropy, c_v and Gibbs and Helmholtz free energy) of the HA system at all qpoints and the whole temperature/pressure range.
Other parameters are the sum of corresponding attributes of all the
Mode objects
. The Helmholtz and Gibbs free energies are defined as:\[ \begin{align}\begin{aligned}F(p,V) = E_{DFT} + F_{vib}(T) = E_{DFT} + U_{vib}(T) - TS(T)\\G(p, V) = F + pV\end{aligned}\end{align} \]- Parameters:
temperature (array[float] | list[float], optional) – Unit: K
pressure (array[float] | list[float], optional) – Unit: GPa
sumphonon (bool) – Whether to sum up the phonon contributions across the sampled q points and take weighted-average.
mutewarning (bool) – Whether print out warning messages of updated temperature and pressure (For QHA).
- Returns:
self.helmholtz (array[float]) – nqpoint*ntemperature. Unit: KJ/mol
self.gibbs (array[float]) – nqpoint*nPressure*nTemperature. Unit: KJ/mol
self.zp_energy (array[float]) – Zero-point energy. nqpoint*1. Unit: KJ/mol
self.u_vib (array[float]) – Vibrational contribution to internal energy. nqpoint*ntemperature. Unit: KJ/mol
self.entropy (array[float]) – nqpoint*ntemperature. Unit: J/mol*K
self.c_v (array[float]) – Constant volume specific heat. nqpoint*ntemperature. Unit: J/mol*K
Note
If
sumphonon = True
, nqpoint = 1.- Raises:
AttributeError – If temperature and pressure are defined neither here nor during initialization
- write_HA_result()
- class Quasi_harmonic(temperature=[], pressure=[], filename=None)
Bases:
object
Generate and rearrange harmonic phonons, store the fitted, volume-dependent QHA phonon information and obtain the QHA thermodynamic properties.
- Parameters:
temperature (array[float] | list[float], optional) – Unit: K
pressure (array[float] | list[float], optional) – Unit: GPa
write_out (bool) – Whether to print the key information into a file.
filename (str) – Name of the output file. Valid if
write_out = True
.
Temperatures and pressures can also be defined by
self.thermodynamics
, whose entries always cover the entries here.Usage:
qha = Quasi_harmonic() qha.from_QHA_file('qha_phonon.out') qha.thermo_freq(eos_method='birch_murnaghan', temperature=[0, 100, 200, 300], pressure=[0., 0.5, 1.]):
- from_HA_files(input_files, scelphono=[], imaginary_tol=-0.0001, q_overlap_tol=0.0001, mode_sort_tol=0.4)
Read data from individual HA calculation outputs. Imaginary modes and overlapped q points are forced to be cleaned.
- Parameters:
input_files (list[str]) – List of phonon output filenames.
scelphono (array[float] | list[float]) – Corresponds to the ‘SCELPHONO’ keyword in CRYSTAL. Either 3*3 or ndimension*ndimension. By default a 1*1*1 ‘SCELPHONO’ is assumed.
imaginary_tol (float) – The threshold of negative frequencies.
q_overlap_tol (float) – The threshold of overlapping points, defined as the 2nd norm of the difference of fractional q vectors
mode_sort_tol (float | None) – The threshold of close mode overlaps. If none, do not sort modes.
- Returns:
self (Quasi_harmonic) – New Attributes listed below
self.ncalc (int) – Number of HA phonon calculations.
self.combined_phonon (list[Harmonic]) – List of Harmonic objects.
self.combined_volume (list[float]) – Volumes. Unit: Angstrom^3
self.combined_edft (list[float]) – DFT total energies. Unit: KJ/mol
self.combined_mode (list[Mode]) – List of mode objects.
- from_QHA_file(input_file, scelphono=[], imaginary_tol=-0.0001, q_overlap_tol=0.0001, mode_sort_tol=0.4)
Read data from a single QHA calculation at Gamma point. Imaginary modes and overlapped q points are forced to be cleaned.
- Parameters:
input_files (str | list[str]) – Only 1 QHA file is permitted.
scelphono (array[float] | list[float])
imaginary_tol (float) – The threshold of negative frequencies.
q_overlap_tol (float) – The threshold of overlapping points, defined as the 2nd norm of the difference of fractional q vectors
mode_sort_tol (float | None) – The threshold of close mode overlaps. If none, do not sort modes.
Returned attributes are consistent with
Quasi_harmonic.from_HA_files
.- Raises:
ValueError – If multiple files are defined.
- from_phonopy_files(phono_yaml, struc_yaml=None, edft=None, imaginary_tol=-0.0001, q_overlap_tol=0.0001, q_id=None, q_coord=None)
Build a QHA object from Phonopy ‘band.yaml’ or ‘qpoints.yaml’ file.
- Parameters:
phono_yaml (list[str]) – ncalc*1 list of Phonopy band.yaml or qpoint.yaml files
struc_yaml (list[str]) – ncalc*1 list of Phonopy phonopy.yaml or phonopy_disp.yaml files. Needed only if a qpoint.yaml file is read.
edft (list[float]) – ncalc*1 list / array of DFT energies.
imaginary_tol (float) – The threshold of negative frequencies.
q_overlap_tol (float) – The threshold of overlapping points, defined as the 2nd norm of the difference of fractional q vectors
q_id (list[int]) – See
Harmonic.from_phonopy
.q_coord (list[list]) – See
Harmonic.from_phonopy
.
Note
q_id
andq_coord
should be set once and are applicable to all the yaml files.Returned attributes are consistent with
Quasi_harmonic.from_HA_files
.
- eos_fit(volume, energy, method, write_out=True, **kwargs)
Fit energy-volume relationship by equation of states.
- Parameters:
volume (array[float]) – Unit: Angstrom^3
energy (array[float]) – Unit: kJ/mol
method (str) – Name of EoS used. Consistent with Pymatgen.
write_out (bool) – Whether to print EOS information.
order (int) – For DeltaFactor / Polynomial methods.
min_ndata_factor (int) – For the NumericalEOS method.
max_poly_order_factor (int) – For the NumericalEOS method.
min_poly_order_factor (int) – For the NumericalEOS method.
- Returns:
eos (Pymatgen EOS) – The fitted equation of state.
eos_method (string) – Name of the fitted equation of state
- freq_polynomial_fit(order)
Fit phonon frequencies as polynomial functions of volumes.
- Parameters:
order (list[int] | array[int]) – The order of polynomials used.
- Returns:
self.fit_order (int) – The optimal order of polynomial fit.
Please also refer to
self.poly_fit
andself.poly_fit_rsquare
attributes of Mode class.
- thermo_freq(eos_method='birch_murnaghan', poly_order=[2, 3], min_method='BFGS', volume_bound=None, mutewarning=False, **kwargs)
Obtain thermodynamic properties by explicitly fitting phonon frequencies as polynomial functions of volume. DFT total energies are fitted as a function of volume by equation of states (EOS).
The equilibrium volume is fitted by minimizing Gibbs free energy at constant temperature and pressure.
\[V(T,p)=\text{min}[G(V;T,p)]=\text{min}[E_{0}(V)+F_{vib}(V;T,p)+pV)]\]- Parameters:
eos_method (str, optional) – EOS used to fit DFT total energy and Helmholtz free energy (to get bulk modules).
poly_order (array[int] | list[int], optional) – The order of polynomials used to fit frequency as the function of volumes.
min_method (string, optional) – Minimisation algorithms.
volume_bound (tuple-like, optional) – equilibrium volumes. Unit: Angstrom^3
mutewarning (bool, optional) – Whether print out warning messages.
temperature (array[float], optional) – Unit: K
pressure (array[float], optional) – Unit: GPa
order (int, optional) – For DeltaFactor / Polynomial EOSs.
min_ndata_factor (int, optional) – For Numerical EOS.
max_poly_order_factor (int, optional) – For Numerical EOS.
min_poly_order_factor (int, optional) – For Numerical EOS.
Note
Valid entries of
eos_method
are consistent with PyMatGen.- Parameterized and tested algorithms for
min_method
: BFGS(no boundary)
L-BFGS-B(with boundary)
- Parameterized and tested algorithms for
- Returns:
self (Quasi_harmonic) – New Attributes listed below
self.temperature (array) – Unit: K
self.pressure (array) – Unit: GPa
self.volume (array) – nPressure*nTemperature, same below. Equilibrium volumes. Unit: \(\AA^{3}\)
self.helmholtz (array) – Helmholtz free energy. Unit: kJ/mol
self.gibbs (array) – Gibbs free energy. Unit: kJ/mol
self.entropy (array) – Entropy. Unit: \(J.mol^{-1}.K^{-1}\)
self.c_v (array) – Constant volume specific heat. Unit: \(J.mol^{-1}.K^{-1}\)
self.e0_eos (Pymatgen EOS) – Pymatgen EOS object. EOS used to fit DFT energy.
self.e0_eos_method (str) – Name of the EOS.
- Raises:
ValueError – If temperature or pressure is defined neither here nor during initialization.
- thermo_gruneisen(eos_method='birch_murnaghan', min_method='BFGS', volume_bound=None, mutewarning=False, **kwargs)
Grüneisen parameters and related properties. The macroscopic Grüneisen parameter is defined as:
\[\gamma=\sum_{\textbf{q}i}\frac{\gamma_{\textbf{q}i}C_{V,\textbf{q}i}}{C_{V}}\]Thermal expansion coefficient in Grüneisen model:
\[\alpha_{V}^{gru}=\frac{\gamma C_{V}}{K_{T}V}\]Note
The Grüneisen model is used to fit frequencies, equivalent to using
self.thermo_freq(poly_order=[1,])
.For arguments, see
self.thermo_freq
.- Returns:
self (Quasi_harmonic) – New attributes listed below. Other attributes are the same as
self.thermo_freq
.self.gruneisen (array) – npressure*ntemperature, same below. Macroscopic Grüneisen parameter. Temperature should > 0.
self.alpha_vgru (array) – Thermal expansion coefficient by Grüneisen method.
self.c_pgru (array) – Constant pressure specific heat by Grüneisen method. Unit: \(J.mol^{-1}.K^{-1}\)
self.k_t (array) – Isothermal bulk modulus. Unit: GPa.
self.k_sgru (array) – Adiabatic bulk modulus by Grüneisen method. Unit: GPa.
- thermo_eos(eos_method='birch_murnaghan', poly_order=[2, 3], mutewarning=False, **kwargs)
Obtain thermodynamic properties by fitting EOS, which is fitted by the Helmholtz free energies of sampled harmonic phonons. The explicit sorting and fitting of frequency-volume relationship is disabled.
Entropy is obtained by taking the derivation of Gibbs free energy at constant pressure.
\[S=-\left(\frac{\partial G}{\partial T}\right)_{p}\]Constant pressure specific heat is obtained by taking the second derivative of \(G\).
\[C_{p}=-T\left(\frac{\partial^{2}G}{\partial T^{2}}\right)_{p}\]Note
poly_order
should >= 2.For arguments, see
self.thermo_freq
.- Returns:
self (Quasi_harmonic) – New attributes listed below
self.temperature (array) – Unit: K
self.pressure (array) – Unit: GPa
self.volume (array) – nPressure*nTemperature, same below. Equilibrium volumes. Unit: \(\AA^{3}\)
self.helmholtz (array) – Helmholtz free energy. Unit: kJ/mol
self.gibbs (array) – Gibbs free energy. Unit: kJ/mol
self.entropy (array) – Entropy. Unit: \(J.mol^{-1}.K^{-1}\)
self.c_p (array) – Constant pressure specific heat. Unit: \(J.mol^{-1}.K^{-1}\)
self.fe_eos (list[Pymatgen EOS]) – nTemperature*1 list of Pymatgen EOS objects. EOSs used to fit HA free energy at constant temperature.
self.fe_eos_method (str) – The name of EOS used.
- Raises:
Exception – If the number of HA calculations is less than 4.
ValueError – If temperature or pressure is defined neither here nor during initialization.
- expansion_vol(poly_order=[2, 3], plot=True, fit_fig='expansion_fit.png')
Fit the thermal expansion curve and get thermal expansion coefficients at equilibrium volumes.
The volumetric thermal expansion coefficient at constant pressure:
\[\alpha_{V}(T) = \frac{1}{V(T)}\left(\frac{\partial V(T)}{\partial T}\right)_{p}\]- Parameters:
poly_order (list[int]) – method = ‘polynomial’, order of polynomials.
plot (bool) – Plot V-T curves to examine the goodness of fitting. An interactive window will pump out to let user to specify the optimial fitting.
fit_fig (str) – File name for fittings. A temperal figure is printed to help the user choose the optimal fitting.
- Returns:
self (Quasi_harmonic) – New attributes listed below
self.vol_fit (list) – nPressure*1. List of Numpy polynomial object, the fitted volume V(T)
self.alpha_v (array) – nPressure*nTemperature. Expansion coefficients at equilibrium volumes
- bulk_modulus(adiabatic=True, **kwargs)
Calculate isothermal and adiabatic bulk moduli at equilibrium volumes.
The following equations are used:
\[ \begin{align}\begin{aligned}K_{T}(p;T) = V(p;T)\left(\frac{\partial^{2}F(V;T)}{\partial V^{2}}\right)_{T}\\K_{S} = K_{T} + \frac{\alpha^{2}_{V}VTK^{2}_{T}}{C_{V}}\end{aligned}\end{align} \]To get \(K_{T}\), Helmholtz free energy is fit as isothermal EOSs. For
self.thermo_eos()
, that means doing nothing; Forself.thermo_freq()
, EOS fitting is required, whose form is the same as EOS used for \(E_{0}\).- Parameters:
adiabatic (bool) – Whether to fit adiabatic bulk modulus. Thermal expansion coefficient needed.
order – (int, optional): To restore EOS.
min_ndata_factor – (int, optional): To restore EOS.
max_poly_order_factor – (int, optional): To restore EOS.
min_poly_order_factor – (int, optional): To restore EOS.
- Returns:
self (Quasi_harmonic) – New attributes listed below
self.k_t (array) – nPressure*nTemperature, same below. Isothermal bulk modulus. Unit: GPa.
self.k_s (array) – Adiabatic bulk modulus. Unit: GPa.
- specific_heat()
Calculate constant volume or pressure specific heat at equilibrium volumes.
The following equation is used:
\[C_{p} - C_{V} = \alpha_{V}^{2}K_{T}VT\]- Returns:
self (Quasi_harmonic) – New attributes listed below
self.c_v (array) – nPressure*nTemperature, same below. Constant volume specific heat. Unit: \(J.mol^{-1}.K^{-1}\)
self.c_p (array) – Constant pressure specific heat. Unit: \(J.mol^{-1}.K^{-1}\)
Note
This method fits
self.c_p
byself.c_v
whenthermo_freq
andthermo_gruneisen
was used.self.c_v
is obtained by whenthermo_eos
is used.
- expansion_lin(poly_order=[2, 3], interp=None)
Fit linear expansions of lattice parameters by the 2-order Taylor expansion.
\[G(\mathbf{p})=G_{0}(\mathbf{p_{0}})+\Delta\mathbf{p}^{T}\mathbf{H}\Delta\mathbf{p}\]:math::G is Gibbs free energy. \(\mathbf{p}\) is the vector of lattice parameters. \(\Delta\mathbf{p}\) means the difference between the fitted and equilibrium lattice parameters. \(\mathbf{H}\) is the Hessian of \(G\) and displacements along lattice parameters.
The RMS deviations (RMSD) of the following equation is minimized at constant temperature and pressure. But deviations from equilibrium volume might occur. RMSD of Gibbs free energy is available in output file only.
\[\mathbf{p_{0}} = \min\left\{\Delta\mathbf{p}^{T}\mathbf{H}\Delta\mathbf{p} - [G(\mathbf{p})-G_{0}(T,p)]\right\}\]This method requires a larger number of HA calculations to ensure a small RMSD. Typically the number of HA calculations should follow the equation below, otherwise the warning massage is given.
\[n_{HA} \geq n_{latt} + \sum_{i=1}^{n_{latt}}i\]\(n_{latt}\) is the lenth of the minimial set of lattice parameters. The optimized lattice parameters at DFT level are used for fitting.
Note
Raimbault, V. Athavale and M. Rossi, Phys. Rev. Materials, 2019, 3, 053605.
- Parameters:
poly_order (list[int]) – Order of polynomials used to fit the linear expansion coefficients. The optimal fit across the sampled temperature and pressure range of a certain lattice parameter is automatically chosen based on \(R^{2}\).
interp (int) – Number of interpolated geometries. All the HA geometries are used besides the interpolated ones.
- Returns:
self (Quasi_harmonic) – New attributes listed below
self.lattice (array) – nPressure*nTemperature*nLattice. The equilibrium values of minimal set of lattice parameters.
self.latt_fit (list) – nPressure*nLattice. Numpy polynomial object, the fitted a(v). Linear part only.
self.alpha_latt (array) – nPressure*nTemperature*nLattice. Linear expansion coefficients. Linear part only.
- class Phonopy
Bases:
object
The convertor between Phonopy and CRYSTALpytools file formats
- classmethod read_structure(file)
Read geometry from Phonopy band.yaml or phonopy.yaml or phonopy_disp.yaml files.
- Parameters:
file (str) – Phonopy yaml file
- Returns:
struc (Pymatgen Structure)
- Raises:
Exception – If the length unit in yaml file is neither ‘au’ nor ‘angstrom’.
- classmethod read_frequency(file, q_id=None, q_coord=None)
Read phonon frequency from Phonopy band.yaml or qpoints.yaml files. Frequency units must be THz (default of Phonopy).
- Parameters:
file (str) – Phonopy yaml file
q_id (list[int]) – Specify the id (from 0) of q points to be read. nqpoint*1 list.
q_coord (list[list]) – Specify the coordinates of q points to be read. nqpoint*3 list.
q_id
andq_coord
should not be set simultaneously. If set,q_id
takes priority andq_coord
is ignored. If both are none, all the points will be read.- Returns:
qpoint (list) – natom*2 list. 1st element: 3*1 array. Fractional coordinates of q points; 2nd element: float. Weight
frequency (array) – nqpint*nmode array. Phonon frequency in THz.
- Raises:
Exception – If (some of) q point is not found.
- classmethod write_force_constants(hessfile='HESSFREQ.DAT', phonopyfile='FORCE_CONSTANTS')
Write Phonopy/VASP FORCE_CONSTANTS file by CRYSTAL HESSFREQ.DAT file.
For example, to convert the calculation ‘example’ with a 4*4*4 supercelland get phonon frequencies at Gamma point, use the following code:
>>> from CRYSTALpytools.thermodynamics import Phonopy >>> Phonopy.write_force_constants(hessfile='example.HESSFREQ') >>> phonopy --crystal --qpoints='0 0 0' -c example.out --dim='4 4 4' --readfc
- Parameters:
hessfile (str) – The HESSFREQ.DAT file
phonopyfile (str) – The output name
- class Output
Bases:
object
Deal with output data file
- classmethod write_HA_result(ha)
Write harmonic phonon information.
- Parameters:
ha (Harmonic) –
CRYSTALpytools.thermodynamic.Harmonic
object
- classmethod write_QHA_combinedata(qha)
Write QHA combined phonon information.
- Parameters:
qha (Quasi_harmonic) –
CRYSTALpytools.thermodynamic.Quasi_harmonic
object
- classmethod write_QHA_sortphonon(qha, close_overlap)
Write QHA phonon sort information.
- Parameters:
qha (Quasi_harmonic) –
CRYSTALpytools.thermodynamic.Quasi_harmonic
object.close_overlap (array[int]) – nqpoint*nmode_ref*ncalc*nmode_sort. Number of close overlaps.
- classmethod write_QHA_eosfit(qha, eos, method)
Write QHA phonon eos fit information.
- Parameters:
qha (Quasi_harmonic) –
CRYSTALpytools.thermodynamic.Quasi_harmonic
object.order (list[int]) – Orders of polynomials
method (str) – Name of EoS used.
- classmethod write_QHA_polyfit(qha, order, rsquare_q)
Write QHA phonon polynomial fit information.
- Parameters:
qha (Quasi_harmonic) –
CRYSTALpytools.thermodynamic.Quasi_harmonic
object.order (list[int]) – List of polynomial orders.
rsquare_q (array) – Nqpoint*Norder list. Overall goodness at a q point.
- classmethod write_QHA_thermofreq(qha, min_method, volume_bound)
Write QHA thermodynamics information (frequency fitting).
- Parameters:
qha (Quasi_harmonic) –
CRYSTALpytools.thermodynamic.Quasi_harmonic
object.
- classmethod write_QHA_thermogru(qha)
Write QHA thermodynamics information (Grüneisen fitting).
- Parameters:
qha (Quasi_harmonic) –
CRYSTALpytools.thermodynamic.Quasi_harmonic
object.
- classmethod write_QHA_thermoeos(qha)
Write QHA thermodynamics information (EOS fitting).
- Parameters:
qha (Quasi_harmonic) –
CRYSTALpytools.thermodynamic.Quasi_harmonic
object.
- classmethod write_expansion_vol(qha, fit_order, fit_rs)
Write volumetric thermal expansions.
- Parameters:
qha (Quasi_harmonic) –
CRYSTALpytools.thermodynamic.Quasi_harmonic
object.fit_order (int) – The order of polynomial used for fitting.
fit_rs (list[float]) – R^2 of fitting. nPressure*1 list.
- classmethod write_bulk_modulus(qha, adiabatic)
Write bulk moduli.
- Parameters:
qha (Quasi_harmonic) –
CRYSTALpytools.thermodynamic.Quasi_harmonic
object.
- adiabatic (bool): Whether the adiabatic bulk modulus \(K_{S}\)
is fitted.
- classmethod write_specific_heat(qha)
Write bulk moduli.
- Parameters:
qha (Quasi_harmonic) –
CRYSTALpytools.thermodynamic.Quasi_harmonic
object.
- classmethod write_expansion_latt(qha, e_err, fit_order, r_square)
Write linear expansion information.
- Parameters:
qha (Quasi_harmonic) –
CRYSTALpytools.thermodynamic.Quasi_harmonic
object.e_err (array) – RMS deviation of Gibbs free energy at T and p.
fit_order (array) – The order of polynomials used for fitting lattice parameter.
r_square (array) – R^2 of fitting. nLatt*nPress
CRYSTALClear.unit_test
- crystal_input_test(folder_path)
- crystal_output_test(folder_path)
- convert_test(folder_path)
- test_all(folder_path)
CRYSTALClear.units
- H_to_eV(energy)
- eV_to_H(energy)
- H_to_kjmol(energy)
- kjmol_to_H(energy)
- au_to_angstrom(length)
- angstrom_to_au(length)
- cm_to_thz(freq)
- thz_to_cm(freq)
- hartree_to_thz(freq)
- thz_to_hartree(freq)
- amu_to_me(mass)
- me_to_amu(mass)