\(\newcommand{\AA}{\text{Å}}\)
CRYSTALpytools.crystal_io module
Objects of input / output files of CRYSTAL. Methods to edit or substract data from corresponding files are provided.
- class Crystal_input(source='')
Bases:
Crystal_inputBASE
Crystal input object inherited from the Crystal_inputBASE object. For the basic set-ups of keywords, please refer to manuals there.
- Parameters:
source (str) – Template file or string. An empty object is generated for
''
.
- classmethod read_file(source)
Instantiate class object from file.
- Parameters:
source (str) – The name of the input file.
- Returns:
cls (Crystal_input)
- write_file(file)
Write data to a file
- 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_user(bs_str, z=[], fmt='crystal', append=False, charge=None, ECP=None, BSfile=None, title='MYBASIS')
A shortcut to get user defined (free format) basis sets from Basis Set Exchange (BSE), formatted string or file.
In some rare cases the user might want to manually change the charge assigned to the atoms to get ionic initial guess.
- Parameters:
bs_str (str) – Name of basis set(BSE), formatted string or file.
z (list[int]) – List of elements, specified by conventional atomic numbers. BSE only.
fmt (str) – Format string. Consistent with BSE python API. For string and files.
append (bool) – Whether to cover old entries. If the old entry contains ‘BASISSET’, it will be removed anyway. Useful when different series of basis sets are defined
charge (dict) – To define charge. Keys: z, conventional atomic number; Values: charge of shells, whose sequence must be consistent with BS definition.
ECP (dict) – Definition of effective core pseudopotentials. Keys: z, conventional atomic number; Values: ECP string.
BSfile (str) – If not None, print basis set definitions into
BSfile
and usetitle
as the entry to ‘BASISSET’ keywordtitle (str) – Useful when
BSfile
is not None.
- bs_keyword(keyword)
A shortcut to define basis sets by keywords. Equivalent to
self.basisset.basisset
- Parameters:
keyword (str) – Name of basis set
- class Crystal_output(output_name=None)
Bases:
object
This class reads a CRYSTAL output and generates an object.
- Parameters:
output_name (str) – Filename
- classmethod read_file(output_name)
Reads a CRYSTAL output file.
- Parameters:
output_name (str) – Name of the output file.
- Returns:
cls (Crystal_output)
- 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(initial=True)
Return the symmetry operators
- Parameters:
initial (bool) – Optimization Only Read symmetry operators of the initial or the final geometry.
- 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 (CStructure | CMolecule) – A modified pymatgen Structure or molecule object.
- 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.
- 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 (CStructure | CMolecule) – A modified 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 wrapper 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
Note
For optimizations, it only returns to optimization convergence. SCF convergence is not available. Call
get_opt_convergence()
if needed.- Parameters:
history (bool) – Deprecated. Always return to history.
- 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_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 in eV.
- get_scf_convergence(all_cycles=False, scflog=None)
Returns the scf convergence energy and energy difference.
Note
With empirical corrections (DFT-D or gCP), the
scf_energy
refers to DFT energy only. Thefinal_energy
attribute includes corrections at the final step.final_energy
is valid only ifscf_status='converged'
.Note
For optimizations, using
get_opt_convergence()
withscf_history=True
is suggested.- Parameters:
all_cycles (bool) – Deprecated.
scflog (str) – Path to ‘SCFOUT.LOG’. For optimizations,
scflog=None
returns to either the initial SCF convergence or every SCF step if ‘ONELOG’ keyword is enabled. This option helps to get SCF history from ‘SCFOUT.LOG’ file.
- Returns:
self (Crystal_output) – New attributes listed below
self.scf_cycles (int|list) – Number of cycles. Array if multiple SCF cycles are read.
self.scf_status (str|list) – ‘terminated’, ‘converged’, ‘too many cycles’ and ‘unknown’. List if multiple SCF cycles are read.
self.scf_energy (array|list) – SCF energy convergence. Unit: eV
self.scf_deltae (array|list) – Energy difference. Unit: eV
self.final_energy (float) – Last step energy with corrections. Unit: eV
- get_fermi_energy(history=False, scflog=None)
Returns the system Fermi energy / valance band maximum (insulators). For
history=True
, Fermi energy at step 0 is always 0 as no Fermi energy is solved. Similary, steps with SPINLOCK or level shifting also return to 0 Fermi energy.- Parameters:
history (bool) – Whether to read the convergence history of Fermi energy.
scflog (str) – Path to ‘SCFOUT.LOG’. For optimizations,
scflog=None
returns to either the initial Fermi energy or energies of every SCF if ‘ONELOG’ keyword is enabled. This option helps to get SCF history from ‘SCFOUT.LOG’ file.
- Returns:
self.fermi_energy (float|array) – Fermi energy of the system. For spin-polarized insulators, the highest VBM is returned.
self.spin_pol (bool) – Not returned but defined Whether the calculation is spin-polarized.
- get_band_gap(history=False, scflog=None)
Returns the system band gap. For
history=True
, gap at step 0 is always 0 as no energy gap is solved. Similary, steps with SPINLOCK or level shifting also return to 0 gap.- Parameters:
history (bool) – Whether to read the convergence history of band gap.
scflog (str) – Path to ‘SCFOUT.LOG’. For optimizations,
scflog=None
returns to either the initial band gap or gaps of every SCF if ‘ONELOG’ keyword is enabled. This option helps to get SCF history from ‘SCFOUT.LOG’ file.
- 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
).self.spin_pol (bool) – Not returned but defined Whether the calculation is spin-polarized.
- 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, \(\alpha\), \(\beta\)].
- get_opt_convergence(primitive=False, scf_history=False, scflog=None, write_gui=False, gui_name=None, symmetry='pymatgen', **kwargs)
Returns optimisation convergence.
- Parameters:
primitive (bool) – Deprecated. Use
get_primitive_geometry()
for initial and last geometries. The method here always returns to the actual structure optimized.scf_history (bool) – Read SCF history of each optimisation step. Also refer to the
get_scf_convergence()
method.scflog (str) – Path to ‘SCFOUT.LOG’. For optimizations,
scflog=None
returns to either the initial SCF convergence or every SCF step if ‘ONELOG’ keyword is enabled. This option helps to get SCF history from ‘SCFOUT.LOG’ file.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) – Modified CStructure/CMolecule 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
self.final_energy (float) – Last step energy with corrections. Unit: eV
- get_forces(initial=True, grad=False, scflog=None)
Read forces.
Note
To get forces of the last optimization step, ‘SCFOUT.LOG’ file or the ‘ONELOG’ keyword are needed.
- Parameters:
initial (bool) – Deprecated.
grad (bool) – Return gradient convergence history. For optimizations only.
scflog (str) – Path to ‘SCFOUT.LOG’. For optimizations,
scflog=None
returns forces from the initial SCF convergence or forces of every SCF step if ‘ONELOG’ keyword is enabled. This option helps to get forces of every step from ‘SCFOUT.LOG’ file.
- 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.
Note
In QHA calculations,
self.nqpoint
refer to harmonic phonons computed. In other cases it refers to actual q points in reciprocal space.Note
This method is developed as a basic I/O function for thermodynamic analysis. For phonon band / dos or IR / Raman spectra, please refer to the
get_phonon_band
,get_phonon_dos
andget_spectra
methods.- Parameters:
read_eigvt (bool) – Whether to read phonon eigenvectors and normalize it to 1.
rm_imaginary (bool) – Set negative frequencies to 0 and remove all the related properties. Only eigenvectors are kept.
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 q points, defined as the 2nd norm of the difference of fractional q vectors
eigvt_amplitude (float|str) – Normalize the eigenvector to a certian amplitude. Either a number or ‘classical’ (classical amplitude in Bohr).
- 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 1*nQpoint list of 1*2 list whose first element is a 1*3 array of fractional coordinates and the second is its weight.
self.nmode (array[int]) – Number of modes at q point. 1*nQpoint array.
self.frequency (array[float]) – nQpoint*nMode array ofvibrational frequency. Unit: THz
self.mode_symm (array[str]) – nQpoint*nMode array of the irreducible representations. In Mulliken symbols.
self.intens (array[float]) – nqpoint*nmode array of harmonic IR 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.
- get_phonon_band(q_overlap_tol=0.0001)
Read phonon bands from CRYSTAL output file.
- Parameters:
q_overlap_tol (float) – The threshold for overlapped k points. Only used for getting tick positions.
- Returns:
self.pband (PhononBand) –
phonons.PhononBand
object.
- get_phonon_dos(read_INS=False, atom_prj=[], element_prj=[])
Read phonon density of states from CRYSTAL output file.
- Parameters:
read_INS (bool) – Read the inelastic neutron scattering spectra.
atom_prj (list) – Read the projections of atoms with specified labels.
element_prj (list) – Read projections of elements with specified conventional atomic numbers.
- Returns:
self.pdos (PhononDOS) –
phonons.PhononDOS
object.
- get_spectra(specfile, type='infer')
Read spectra from IRSPEC.DAT / RAMSPEC.DAT files.
Note
In principle, the code cannot tell the type of the files. It is important to assign the appropriate type to the
type
keyword. Currently the available options are ‘IRSPEC’ and ‘RAMSPEC’.- Parameters:
specfile (str) – File name of spectra data.
type (str) – Type of the file. See above. If ‘infer’, type is inferred from input file name.
- Returns:
self.* (IR|Raman) – Dependending on the
type
keyword, return tospectra.IR
orspectra.Raman
objects. Attribute names same astype
are set.
- get_elatensor(*thickness)
Extracts the elastic tensor from the data.
- Parameters:
*thickness (float) – Effective thickness of low dimensional materials, in \(\AA\).
- Returns:
self.tensor (numpy.ndarray) – Symmetrized elastic tensor in Voigt notation. For 3D systems, 6*6; for 2D, 3*3; for 1D, 1*1. The matrix always has 2 dimensions. Unit: GPa for 3D and 1D, 2D if effective thickness of materials are specified. Otherwise GPa.m for 2D and GPa.m \(^{2}\) for 1D (might lead to very small numbers).
- get_anh_spectra()
Note
This is not for the released feature of CRYSTAL23 v1.0.1
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.
- property n_atoms
Deprecated. Get structure object by ‘get_geometry’ and call the ‘num_sites’ attribute
- property atom_symbols
Deprecated. Get structure object by ‘get_geometry’ and call the ‘species_symbol’ attribute
- property atom_numbers
Deprecated. Get structure object by ‘get_geometry’ and call the ‘species_Z’ attribute
- property atom_positions
Deprecated. Get structure object by ‘get_geometry’ and call the ‘crys_coords’ attribute
- property atom_positions_cart
Deprecated. Get structure object by ‘get_geometry’ and call the ‘cart_coords’ attribute
- property atom_positions_frac
Deprecated. Get structure object by ‘get_geometry’ and call the ‘crys_coords’ attribute
- get_last_geom(write_gui_file=True, symm_info='pymatgen')
Deprecated. Use
get_geometry(initial=False)
.
- class Properties_input(source='')
Bases:
Properties_inputBASE
Properties input object inherited from the Properties_inputBASE object. For the basic set-ups of keywords, please refer to manuals there.
- Parameters:
source (str) – Template file or string. An empty object is generated for
''
.
- classmethod read_file(source)
Instantiate the object from a file.
- Parameters:
source (str) – The name of the input file.
- Returns:
cls (Properties_input)
- write_file(file)
Writes the properties input to a file.
- Parameters:
file (str) – The name of the d3 file.
- make_band_block(k_path, n_kpoints, first_band, last_band, print_eig=0, print_option=1, precision=5, title='BAND STRUCTURE CALCULATION')
A shortcut for
self.band()
. Also accepts pymatgen HighSymmKpath objects for k path.- Parameters:
k_path (list | 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) – Precision in the calculation of
numpy.gcd
title (str) – The title of the calculation.
- make_doss_block(n_points=200, band_range=None, e_range=None, plotting_option=2, poly=12, print_option=1, projections=[], proj_type='atom', output_file=None)
A shortcut for
self.doss()
.- Parameters:
n_points (int) – The number of points in the DOS plot.
band_range (list or tuple) – The range of bands to include in the DOS calculation.
e_range (list or tuple) – The energy range for the DOS calculation. Unit: eV
plotting_option (int) – DOS plotting options.
poly (int) – Degree of the polynomial for smearing.
print_option (int) – Properties printing options.
projections (list) – nProj*nElement list of ints specifying the projections for the pdoss calculation, or nProj*1 list of strings specifying element names to project on all atoms of the same element (
output_file
is needed).proj_type (str) – Type of projection (‘atom’ or ‘site’).
output_file (str) – Output file of ‘crystal’ calculation.
- 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)
Deprecated. Use either
self.doss()
orself.make_doss_block()
- make_newk_block(shrink1, shrink2, Fermi=1, print_option=0)
Deprecated. Use
self.newk()
.
- make_bands_block(k_path, n_kpoints, first_band, last_band, print_eig=0, print_option=1, precision=5, title='BAND STRUCTURE CALCULATION')
Deprecated. Use
self.make_band_block()
.
- from_file(input_name)
Deprecated. Use
read_file()
- write_properties_input(input_name)
Deprecated. Use
write_file()
- class Properties_output(properties_output=None)
Bases:
POutBASE
Creates a Properties_output object. Since not all results from ‘properties’ executable creates are summarized in output, methods of this class usually requires multiple files.
Note
Most of methods directly dealing with output file itself are saved in
base.output.POutBASE
object. Users typically do not need to call methods there since they are called internally.- Parameters:
properties_output (str) – Properties output file
- classmethod read_file(properties_output)
Parse the properties output file.
- Parameters:
properties_output (str) – The properties output file.
- Returns:
cls (Properties_output)
- read_electron_band(band_file)
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
- Returns:
self.bands (ElectronBand) – An
CRYSTALpytools.electronics.ElectronBand
object
- read_electron_dos(dos_file)
Generate doss object from CRYSTAL DOSS.DAT or fort.25 file. Energy unit: eV. E Fermi is aligned to 0. :param dos_file: Name of DOSS.DAT or fort.25 file :type dos_file: str
- Returns:
self.doss (ElectronDOS) – An
CRYSTALpytools.electronics.ElectronDOS
object
- read_topond(topondfile, type='infer')
Read the 2D scalar plot files (‘SURF*.DAT’) or trajectory files (TRAJ*.DAT) written by TOPOND.
Geometry information is printed in the standard ouput, which is not mandatory for ‘SURF*.DAT’ but is mandatory for ‘TRAJ*.DAT’
Note
For the convenience of analysis and plotting, it is important to select the correct type for your input file. By default type=’infer’ will search for (case insensitive) the following strings:
‘SURFRHOO’, ‘SURFSPDE’, ‘SURFLAPP’, ‘SURFLAPM’, ‘SURFGRHO’, ‘SURFKKIN’, ‘SURFGKIN’, ‘SURFVIRI’, ‘SURFELFB’, ‘TRAJGRAD’, ‘TRAJMOLG’.
For their meanings, please refer TOPOND manual.
- Parameters:
topondfile (str) – TOPOND formatted 2D plot file
type (str) – ‘infer’ or specified. Otherwise warning will be given.
- Returns:
self.* (ChargeDensity|SpinDensity|Gradient|Laplacian|HamiltonianKE|LagrangianKE|VirialField|ELF|GradientTraj|ChemicalGraph) – Return to
topond
property classes, depending on input file types. The attribute name is upper case type names. If unknown, return toself.TOPOND
. Atopond.ChargeDensity
ortopond.GradientTraj
class is generated.
- read_ECHG(*f25_files, method='normal', index=None)
Read charge / spin density data from a file. Unit: \(e.\AA^{-3}\).
Available methods are:
‘normal’: Normal 1-file reading.
- ‘substact’: Substracting data from the first entry based on following
entries. Multiple entries or 1 entry with ‘PATO’ keyword enabled. For multiple entries, make sure the charge map is in the first (and ideally the only) ‘MAPN’ data block, otherwise the code won’t get the correct data. For 1 entry with ‘PATO’, data will be substracted from the ‘normal’ system.
- ‘alpha_beta’: Save spin-polarized data in \(\alpha\) /
\(\beta\) states, rather than charge(\(\alpha+\beta\)) / spin(\(\alpha-\beta\)). Single entry only.
Note
The standard screen output is highly recommended to add, which identifies the indices of corresponding 2D data maps. Otherwise the
index
input can be specified. Otherwise, the code only reads the first 2D data map for spin =1 and first 2 maps for spin=2.- Parameters:
*f25_files (str) – Path to the fort.25 file(s).
method (str) – Data processing method. See above.
index (int|list) – Sequence number of headers with the ‘-%-MAPN’ pattern. Useful only if the standard screen output is not available. Starting from 0.
- Returns:
self.echg (ChargeDensity) –
electronics.ChargeDensity
object.
- read_relativistics(f25_file, type, index=None)
Read 2D scalar / vector fields from 2c-SCF calculations, generated by the ‘PROPS2COMP’ keyword.
Note
The standard screen output is highly recommended to add, which identifies the indices of corresponding 2D data maps. Otherwise the
index
must be specified by integer or list of integers.- Parameters:
f25_file (str) – File name of the fort.25 file.
type (str|int|list) – Property to calculate, ‘DENSITY’, ‘MAGNETIZ’, ‘ORBCURDENS’, or ‘SPICURDENS’.
index (int|list) – Sequence number of headers with the ‘-%-MAPN’ pattern. Useful only if the standard screen output is not available. Starting from 0.
- Returns:
self.* (ChargeDensity|Magnetization|OrbitalCurrentDensity|SpinCurrentDensity) – Return to classes defined in the
relativisitcs
module, depending ontype
. The attribute name is upper case type names. Unit: charge densities, \(\AA^{-3}\); magnetization, A/m; Orbital/spin densities, A/m \(^{2}\).
- read_XRDspec(option='LP')
Read XRD spectra from standard screen output of properties calculation. It is envoked by keyword ‘XRDSPEC’.
- Parameters:
option (str) – ‘NC’ for no correction (The ‘INTENS’ col); ‘LP’ for Lorentz and polarization effects (‘INTENS-LP’) and ‘DW’ for LP with Debye-Waller thermal factors (‘INTENS-LP-DW’).
- Returns:
self.XRDspec (XRD) – The
spectra.XRD
object with spectra information.
- 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_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_transport(boltztra_out)
Read electron transport properties by the BOLTZTRA keyword, including ‘KAPPA’, ‘SIGMA’, ‘SIGMAS’, ‘SEEBECK’ and ‘TDF’. Though currently the geometry information is not required, it is saved if the standard output file is given.
Note
For ‘SEEBECK’, all the 9 elements of the tensor was printed. As far as the developers have been aware of, it is symmetrized. Therefore the redundant ‘yx’, ‘zx’ and ‘zy’ dimensions are removed to keep consistent with other outputs.
- Parameters:
boltztra_out (str) – ‘DAT’ files by CRYSTAL BOLTZTRA keyword.
- Returns:
self.* (Tensor|Distribution) –
transport.Tensor
(‘KAPPA’, ‘SIGMA’, ‘SIGMAS’, ‘SEEBECK’) ortransport.Distribution
(TDF) classes, depending on the input file. The attribute name is upper case types.
- read_cry_band(band_file)
Deprecated. Use
read_electron_band
.
- read_cry_doss(dos_file)
Deprecated. Use
read_electron_dos
.
- read_cry_ECHG(f25_file)
Deprecated. Use
read_ECHG
.
- read_cry_ECHG_delta(f25_file1, f25_file2)
Deprecated. Use
read_ECHG
.
- read_cry_contour(properties_output)
Deprecated. Use
read_topond
.
- read_cry_seebeck(properties_output)
Deprecated. Use
read_transport
.
- read_cry_sigma(properties_output)
Deprecated. Use
read_transport
.
- read_vecfield(properties_output, which_prop)
Deprecated. Use
read_relativistics
.
- read_cry_xrd_spec(properties_output)
Deprecated. Use
read_XRDspec
.
- class Crystal_gui(dimensionality=None, lattice=None, symmops=None, atom_number=None, atom_positions=None, space_group=None)
Bases:
Crystal_gui
Inherited from
geometry.Crystal_gui
. See documentations there.
- 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.