Quantum ESPRESSO
Tools for the quantum espresso software package (adaptor). It provide a lot of functions to deal with the qe-related files.
Source: qe.jl
Contents
Index
ZenCore.AtomicPosition
ZenCore.AtomicPositionsCard
ZenCore.AtomicSpecies
ZenCore.AtomicSpeciesCard
ZenCore.AutoKmeshCard
ZenCore.GammaPointCard
ZenCore.KPointsCard
ZenCore.MonkhorstPackGrid
ZenCore.QECard
ZenCore.QEInputEntry
ZenCore.QENamelist
ZenCore.ReciprocalPoint
ZenCore.SpecialPointsCard
Base.delete!
Base.getindex
Base.haskey
Base.parse
Base.setindex!
Base.tryparse
Base.write
ZenCore.adaptor_call
ZenCore.dft_call
ZenCore.dft_resume
ZenCore.dft_stop
ZenCore.qe_adaptor
ZenCore.qe_exec
ZenCore.qe_init
ZenCore.qe_save
ZenCore.qe_to_plo
ZenCore.qe_to_wan
ZenCore.qec_input
ZenCore.qeio_band
ZenCore.qeio_eigen
ZenCore.qeio_energy
ZenCore.qeio_fermi
ZenCore.qeio_kmesh
ZenCore.qeio_lattice
ZenCore.qeq_files
Functions
ZenCore.dft_call
— Functiondft_call(::VASPEngine, it::IterInfo)
Try to carry out full DFT calculation with the vasp
code. It is only a dispatcher. Similar function is defined in qe.jl
as well.
See also: _engine_
.
dft_call(::QEEngine, it::IterInfo)
Try to carry out full DFT calculation with the qe
code. It is only a dispatcher. Similar function is defined in vasp.jl
as well.
See also: _engine_
.
ZenCore.dft_stop
— Functiondft_stop(::VASPEngine)
Try to terminate DFT calculation and kill running process of the DFT backend. It supports the vasp
code. It is only a dispatcher. Similar function is defined in qe.jl
as well.
See also: _engine_
.
dft_stop(::QEEngine)
Try to terminate DFT calculation and kill running process of the DFT backend. It supports the qe
code. It is only a dispatcher. Similar function is defined in vasp.jl
as well.
See also: _engine_
.
ZenCore.dft_resume
— Functiondft_resume(::VASPEngine)
Try to wake up the DFT backend and resume the DFT calculation. It only supports the vasp
code. It is only a dispatcher. Similar function is defined in qe.jl
as well.
See also: _engine_
.
dft_resume(::QEEngine)
Try to wake up the DFT backend and resume the DFT calculation. It only supports the qe
code. It is only a dispatcher. Similar function is defined in vasp.jl
as well.
See also: _engine_
.
ZenCore.adaptor_call
— Functionadaptor_call(::VASPEngine, D::Dict{Symbol,Any})
It is a dispatcher for the DFT-DMFT adaptor. It calls vasp_adaptor()
function to deal with the outputs of the DFT backend (such as vasp) and generate key dataset for the next level adaptor (PLOAdaptor
). Note that similar function is also defined in qe.jl
.
adaptor_call(::QEEngine, D::Dict{Symbol,Any})
It is a dispatcher for the DFT-DMFT adaptor. It calls qe_adaptor()
function to deal with the outputs of the DFT backend (such as qe) and generate key dataset for the next level adaptor (WANNIERAdaptor
). Note that similar function is also defined in vasp.jl
.
See also: qe_adaptor
.
adaptor_call(::PLOAdaptor,
D::Dict{Symbol,Any},
ai::Array{Impurity,1})
It is a dispatcher for the DFT-DMFT adaptor. It calls plo_adaptor()
function to deal with the outputs of the DFT backend (such as vasp) and generate key dataset for the IR adaptor. Note that similar functions are also defined in vasp.jl
, qe.jl
, and wannier.jl
.
See also: plo_adaptor
.
adaptor_call(::WANNIERAdaptor,
D::Dict{Symbol,Any},
ai::Array{Impurity,1})
It is a dispatcher for the DFT-DMFT adaptor. It calls wannier_adaptor()
function to deal with the outputs of wannier90
code and generate key dataset for the IR adaptor. Note that similar functions are also defined in vasp.jl
, qe.jl
, and plo.jl
.
See also: wannier_adaptor
.
ZenCore.qe_adaptor
— Functionqe_adaptor(D::Dict{Symbol,Any})
Adaptor support for the quantum espresso
(pwscf
code). It will parse the output files of the quantum espresso
code, extract the Kohn-Sham dataset, and then fulfill the DFTData
dict (i.e D
).
The following output files of the quantum espresso
code are needed:
scf.out
nscf.out
Note in the input file of the quantum espresso
code, the verbosity
parameter must be set to 'high'.
See also: wannier_adaptor
, ir_adaptor
.
ZenCore.qe_to_wan
— Functionqe_to_wan(D::Dict{Symbol,Any})
Try to call the wannier90
and pw2wannier90
codes to generate the maximally-localized wannier functions. the DFTData
dict (i.e D
) will not be modified.
See also: wannier_adaptor
.
ZenCore.qe_to_plo
— Functionqe_to_plo(D::Dict{Symbol,Any})
Postprocess outputs of the quantum espresso
(pwscf
code), call the wannier90
and pw2wannier90
codes to generate the projected local orbitals (which are not maximally-localized wannier functions). The key data are fed into the DFTData
dict (i.e D
).
Most of the functions used in the qe_to_plo()
function are implemented in another file (wannier.jl
).
See also: plo_adaptor
, qe_adaptor
.
ZenCore.qe_init
— FunctionZenCore.qe_exec
— Functionqe_exec(it::IterInfo, scf::Bool = true)
Execute the quantum espresso
(pwscf
) program, monitor the convergence progress, and output the relevant information. The argument scf
controls which input file should be used. If scf == true
, then the input file is case.scf
, or else it is case.nscf
.
In order to execute this function correctly, you have to setup the following environment variables:
- QE_HOME
and make sure the file MPI.toml
is available.
ZenCore.qe_save
— FunctionZenCore.qec_input
— Functionqec_input(it::IterInfo)
It will parse the QE.INP
file at first. The QE.INP
is a standard, but mini input file for quantum espresso
(pwscf
). It only includes three namelists (namely control
, system
, and electrons
) and three cards (namely ATOMIC_SPECIES
, ATOMIC_POSITIONS
, and K_POINTS
). If you want to support more input entries, please make your own modifications here.
Then this function will try to customize these namelists and cards according to the setup in case.toml
.
At last, it will try to generate the input files for quantum espresso
(pwscf
). They are case.scf
and case.nscf
. As shown by their names, one file is for the self-consistent calculation, while another one is for the non-self-consistent calculation.
The return values of this function are namelist (control
) and card (ATOMIC_SPECIES
), which will be used to check the pseudopotentials within the qe_init()
function.
See also: QENamelist
, QECard
.
ZenCore.qeq_files
— Functionqeq_files(f::String)
Check the essential output files by quantum espresso
(pwscf
). Here f
means only the directory that contains the desired files.
See also: adaptor_core
.
qeq_files()
Check the essential output files by quantum espresso
(pwscf
) in the current directory.
See also: adaptor_core
.
ZenCore.qeio_energy
— Functionqeio_energy(f::String)
Reading quantum espresso's scf.out
file, return DFT total energy, which will be used to determine the DFT + DMFT energy. Here f
means only the directory that contains scf.out
.
qeio_energy()
Reading quantum espresso's scf.out
file, return DFT total energy, which will be used to determine the DFT + DMFT energy.
ZenCore.qeio_lattice
— Functionqeio_lattice(f::String, silent::Bool = true)
Reading quantum espresso's scf.out
file, return crystallography information. Here f
means only the directory that contains scf.out
.
See also: Lattice
, irio_lattice
.
qeio_lattice()
Reading quantum espresso's scf.out
file, return crystallography information.
See also: Lattice
, irio_lattice
.
ZenCore.qeio_kmesh
— Functionqeio_kmesh(f::String)
Reading quantum espresso's nscf.out
file, return kmesh
and weight
. Here f
means only the directory that contains nscf.out
.
Note in scf.out
, the 𝑘-mesh is not uniform. So we have to read k-mesh from the nscf.out
. In addition, the verbosity parameter must be set to 'high' in the input file.
See also: irio_kmesh
.
ZenCore.qeio_eigen
— Functionqeio_eigen(f::String)
Reading quantum espresso's nscf.out
file, return energy band structure information. Here f
means only the directory that contains nscf.out
.
Note that in scf.out
, the eigenvalues may be not defined on the uniform 𝑘-mesh. So we have to read eigenvalues from the nscf.out
file.
Note that the eigenvalues read from nscf.out
is somewhat coarse. They should be updated by the values read from case.eig
.
See also: irio_eigen
.
qeio_eigen()
Reading quantum espresso's nscf.out
file, return energy band structure information.
See also: irio_eigen
.
ZenCore.qeio_fermi
— Functionqeio_fermi(f::String, silent::Bool = true)
Reading quantum espresso's nscf.out
file, return the fermi level. Here f
means only the directory that contains scf.out
.
See also: irio_fermi
.
ZenCore.qeio_band
— Functionqeio_band(f::String)
Reading quantum espresso's bands.out
file, return energy band structure information. Here f
means only the directory that contains bands.out
.
In order to generate the bands.out
file, please change the calculation mode of quantum espresso to bands
, and redirect the standard output to the bands.out
file.
Note that the difference between qeio_band()
and qeio_eigen()
is that the former does not return the occupy matrix. This function should not be called in the DFT + DMFT iterations.
qeio_band()
Reading quantum espresso's bands.out
file, return energy band structure information.
ZenCore.ReciprocalPoint
— TypeReciprocalPoint
Represent a special point of the 3D Brillouin zone. Each of them has a weight w
.
Members
- coord -> Coordinates, i.e., $k_x$, $k_y$, and $k_z$.
- weight -> Weight for the $k$-point.
See also: MonkhorstPackGrid
.
ZenCore.MonkhorstPackGrid
— TypeMonkhorstPackGrid
Represent the Monkhorst-Pack grid.
Members
- mesh -> A length-three vector specifying the $k$-point grid ($nk_1 × nk_2 × nk_3$) as in Monkhorst-Pack grids.
- shift -> A length-three vector specifying whether the grid is displaced by half a grid step in the corresponding directions.
See also: ReciprocalPoint
.
ZenCore.AtomicSpecies
— TypeAtomicSpecies
Represent each line of the ATOMIC_SPECIES
card in the input file of quantum espresso
(pwscf
).
Members
- atom -> Label of the atom. Maximum total length cannot exceed 3 characters.
- mass -> Mass of the atomic species in atomic unit. Used only when performing molecular dynamics (MD) run or structural optimization runs using damped MD.
- upf -> File containing pseudopotential for this species.
Examples
julia> AtomicSpecies("C1", 12, "C.pbe-n-kjpaw_psl.1.0.0.UPF")
AtomicSpecies("C1", 12.0, "C.pbe-n-kjpaw_psl.1.0.0.UPF")
See also: AtomicSpeciesCard
.
ZenCore.AtomicPosition
— TypeAtomicPosition
Represent each line of the ATOMIC_POSITIONS
card in the input file of quantum espresso
(pwscf
).
Members
- atom -> Label of the atom as specified in
AtomicSpecies
. It accepts at most 3 characters. - pos -> Atomic positions. A three-element vector of floats.
- ifpos -> Component
i
of the force for this atom is multiplied by `ifpos(i), which must be either
0or
1`. Used to keep selected atoms and/or selected components fixed in MD dynamics or structural optimization run.
Examples
julia> AtomicPosition('O', [0, 0, 0])
AtomicPosition("O", [0.0, 0.0, 0.0], Bool[1, 1, 1])
See also: AtomicPositionsCard
.
ZenCore.QEInputEntry
— TypeQEInputEntry
An abstract type representing an input component of quantum espresso
(pwscf
). Note that all other input types (such as QENamelist
and QECard
) should subtype QEInputEntry
. It is used to build a internal type system.
See also: QECard
, QENamelist
.
ZenCore.QENamelist
— TypeQENamelist
Represent a namelist in the input file of quantum espresso
(pwscf
), a basic Fortran data structure.
Members
- name -> Name of the namelist. It should be
control
,system
, orelectrons
. If you want to support more namelists, please make your own modifications. - data -> A dict containing pairs of key and value.
See also: QECard
.
ZenCore.QECard
— TypeQECard
It represents abstract cards in the input file of quantum espresso
(pwscf
). It is used to build the internal type system. The input file of quantum espresso
(pwscf
) consists of various cards and namelists, represented by QECard
and QENamelist
, respectively.
See also: QENamelist
.
ZenCore.KPointsCard
— TypeKPointsCard
Represent an abstract card (K-POINTS
) in the input file of quantum espresso
(pwscf
).
ZenCore.AtomicSpeciesCard
— TypeAtomicSpeciesCard
Represent the ATOMIC_SPECIES
card in the input file of quantum espresso
(pwscf
).
Members
- data -> A vector containing
AtomicSpecies
.
See also: AtomicSpecies
.
ZenCore.AtomicPositionsCard
— TypeAtomicPositionsCard
Represent the ATOMIC_POSITIONS
card in the input file of quantum espresso
(pwscf
).
Members
- data -> A vector containing
AtomicPosition
s. - option -> The scheme about how to define atomic positions.
See also: AtomicPosition
.
ZenCore.AutoKmeshCard
— TypeAutoKmeshCard
Represent the K_POINTS
card in the input file of quantum espresso
(pwscf
) (be compatible with the automatic
mode only).
See also: KPointsCard
.
ZenCore.GammaPointCard
— TypeGammaPointCard
Represent the K_POINTS
card in the input file of quantum espresso
(pwscf
) (be compatible with the gamma
mode).
See also: KPointsCard
.
ZenCore.SpecialPointsCard
— TypeSpecialPointsCard
Represent the K_POINTS
card in the input file of quantum espresso
(pwscf
) (be compatible with the tpiba
or crystal
mode).
Members
- data -> A vector containing
ReciprocalPoint
s. - option -> The way about how to define $k$-mesh.
See also: KPointsCard
.
Base.haskey
— FunctionBase.haskey(qnl::QENamelist, key::AbstractString)
Examine the existence of an entry (specified by key
) in the namelist object (qnl
).
See also: QENamelist
.
Base.getindex
— FunctionBase.getindex(qnl::QENamelist, key::AbstractString)
Return an entry (specified by key
) in the namelist object (qnl
).
See also: QENamelist
.
Base.setindex!
— FunctionBase.setindex!(qnl::QENamelist, value, key::AbstractString)
Modify an entry (specified by key
) in the namelist object (qnl
).
See also: QENamelist
.
Base.delete!
— FunctionBase.delete!(qnl::QENamelist, key::AbstractString)
Remove an entry (specified by key
) in the namelist object (qnl
).
See also: QENamelist
.
Base.tryparse
— FunctionBase.tryparse(::Type{AtomicSpeciesCard}, str::AbstractString)
Try to parse the AtomicSpeciesCard
object.
See also: AtomicSpeciesCard
.
Base.tryparse(::Type{AtomicPositionsCard}, str::AbstractString)
Try to parse the AtomicPositionsCard
object.
See also: AtomicPositionsCard
.
Base.tryparse(::Type{KPointsCard}, str::AbstractString)
Try to parse the KPointsCard
object.
See also: KPointsCard
.
Base.tryparse(::Type{AutoKmeshCard}, str::AbstractString)
Try to parse the AutoKmeshCard
object.
See also: AutoKmeshCard
.
Base.tryparse(::Type{GammaPointCard}, str::AbstractString)
Try to parse the GammaPointCard
object.
See also: GammaPointCard
.
Base.tryparse(::Type{SpecialPointsCard}, str::AbstractString)
Try to parse the SpecialPointsCard
object.
See also: SpecialPointsCard
.
Base.parse
— FunctionBase.parse(::Type{QENamelist}, strs::Vector{String}, name::String)
Parse the QENamelist
object. The name of the namelist is specified by argument name
.
See also: QENamelist
.
Base.parse(::Type{T}, str::AbstractString)
Parse the QECard
object. Now we support the following cards:
ATOMIC_SPECIES
(AtomicSpeciesCard
)ATOMIC_POSITIONS
(AtomicPositionsCard
)K_POINTS
(AutoKmeshCard
,GammaPointCard
,SpecialPointsCard
)
See also: QECard
.
Base.write
— FunctionBase.write(io::IO, x::AtomicSpeciesCard)
Write the AtomicSpeciesCard
object to IOStream
.
See also: AtomicSpeciesCard
.
Base.write(io::IO, x::AtomicPositionsCard)
Write the AtomicPositionsCard
object to IOStream
.
See also: AtomicPositionsCard
.
Base.write(io::IO, x::AutoKmeshCard)
Write the AutoKmeshCard
object to IOStream
.
See also: AutoKmeshCard
.
Base.write(io::IO, x::GammaPointCard)
Write the GammaPointCard
object to IOStream
.
See also: GammaPointCard
.
Base.write(io::IO, x::SpecialPointsCard)
Write the SpecialPointsCard
object to IOStream
.
See also: SpecialPointsCard
.