QuantumESPRESSOBase.PWscf module
Types
CrystallographyBase.Lattice — TypeLattice(nml::SystemNamelist)Create a Lattice from a SystemNamelist.
Lattice(card::CellParametersCard)Create a Lattice from a CellParametersCard.
Lattice(card::PWInput)Create a Lattice from a PWInput.
QuantumESPRESSOBase.PWscf.ControlNamelist — TypeControlNamelist(calculation, title, verbosity, restart_mode, wf_collect, nstep, iprint, tstress, tprnfor, dt, outdir, wfcdir, prefix, lkpoint_dir, max_seconds, etot_conv_thr, forc_conv_thr, disk_io, pseudo_dir, tefield, dipfield, lelfield, nberrycyc, lorbm, lberry, gdir, nppstr, lfcpopt, gate)
ControlNamelist(; kwargs...)
ControlNamelist(::ControlNamelist; kwargs...)Represent the CONTROL namelist of pw.x.
QuantumESPRESSOBase.PWscf.SystemNamelist — TypeSystemNamelist(ibrav, celldm, A, B, C, cosAB, cosAC, cosBC, nat, ntyp, nbnd, tot_charge, starting_charge, tot_magnetization, starting_magnetization, ecutwfc, ecutrho, ecutfock, nr1, nr2, nr3, nr1s, nr2s, nr3s, nosym, nosym_evc, noinv, no_t_rev, force_symmorphic, use_all_frac, occupations, one_atom_occupations, starting_spin_angle, degauss, smearing, nspin, noncolin, ecfixed, qcutz, q2sigma, input_dft, exx_fraction, screening_parameter, exxdiv_treatment, x_gamma_extrapolation, ecutvcut, nqx1, nqx2, nqx3, localization_thr, lda_plus_u, lda_plus_u_kind, Hubbard_U, Hubbard_J0, Hubbard_alpha, Hubbard_beta, starting_ns_eigenvalue, U_projection_type, edir, emaxpos, eopreg, eamp, angle1, angle2, constrained_magnetization, fixed_magnetization, lambda, report, lspinorb, assume_isolated, esm_bc, esm_w, esm_efield, esm_nfit, fcp_mu, vdw_corr, london, london_s6, london_c6, london_rvdw, london_rcut, ts_vdw_econv_thr, ts_vdw_isolated, xdm, xdm_a1, xdm_a2, space_group, uniqueb, origin_choice, rhombohedral, zgate, relaxz, block, block_1, block_2, block_height)
SystemNamelist(; kwargs...)
SystemNamelist(::SystemNamelist; kwargs...)Represent the SYSTEM namelist of pw.x.
QuantumESPRESSOBase.PWscf.ElectronsNamelist — TypeElectronsNamelist(electron_maxstep, scf_must_converge, conv_thr, adaptive_thr, conv_thr_init, conv_thr_multi, mixing_mode, mixing_beta, mixing_ndim, mixing_fixed_ns, diagonalization, ortho_para, diago_thr_init, diago_cg_maxiter, diago_david_ndim, diago_full_acc, efield, efield_cart, efield_phase, startingpot, startingwfc, tqr)
ElectronsNamelist(; kwargs...)
ElectronsNamelist(::ElectronsNamelist; kwargs...)Represent the ELECTRONS namelist of pw.x.
QuantumESPRESSOBase.PWscf.IonsNamelist — TypeIonsNamelist(ion_dynamics, ion_positions, pot_extrapolation, wfc_extrapolation, remove_rigid_rot, ion_temperature, tempw, tolp, delta_t, nraise, refold_pos, upscale, bfgs_ndim, trust_radius_max, trust_radius_min, trust_radius_ini, w_1, w_2)
IonsNamelist(; kwargs...)
IonsNamelist(::IonsNamelist; kwargs...)Represent the IONS namelist of pw.x.
Input this namelist only if calculation is "relax", "md", "vc-relax", or "vc-md".
QuantumESPRESSOBase.PWscf.CellNamelist — TypeCellNamelist(cell_dynamics, press, wmass, cell_factor, press_conv_thr, cell_dofree)
CellNamelist(; kwargs...)
CellNamelist(::CellNamelist; kwargs...)Represent the CELL namelist of pw.x.
Input this namelist only if calculation is "vc-relax" or "vc-md".
QuantumESPRESSOBase.PWscf.DosNamelist — TypeDosNamelist(prefix, outdir, ngauss, degauss, Emin, Emax, DeltaE, fildos)
DosNamelist(; kwargs...)
DosNamelist(::DosNamelist; kwargs...)Represent the DOS namelist of dos.x.
QuantumESPRESSOBase.PWscf.BandsNamelist — TypeBandsNamelist(prefix, outdir, filband, spin_component, lsigma, lp, filp, lsym, no_overlap, plot_2d, firstk, lastk)
BandsNamelist(; kwargs...)
BandsNamelist(::BandsNamelist; kwargs...)Represent the BANDS namelist of bands.x.
QuantumESPRESSOBase.PWscf.AtomicSpecies — TypeAtomicSpecies(atom::Union{Char,String}, mass, pseudopot)
AtomicSpecies(x::AtomicPosition, mass, pseudopot)Represent each line in the ATOMIC_SPECIES card in pw.x input files.
The atom field accepts no more than 3 characters.
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")
julia> AtomicSpecies(
AtomicPosition('S', [0.500000000, 0.288675130, 1.974192764]),
32.066,
"S.pz-n-rrkjus_psl.0.1.UPF",
)
AtomicSpecies("S", 32.066, "S.pz-n-rrkjus_psl.0.1.UPF")QuantumESPRESSOBase.PWscf.AtomicSpeciesCard — TypeAtomicSpeciesCard <: CardRepresent the ATOMIC_SPECIES card in QE. It does not have an "option".
QuantumESPRESSOBase.PWscf.AtomicPosition — TypeAtomicPosition(atom::Union{Char,String}, pos[, if_pos])
AtomicPosition(x::AtomicSpecies, pos, if_pos)Represent each line in the ATOMIC_POSITIONS card in pw.x input files.
The atom field accepts no more than 3 characters.
Examples
julia> AtomicPosition('O', [0, 0, 0])
AtomicPosition("O", [0.0, 0.0, 0.0], Bool[1, 1, 1])
julia> AtomicPosition(
AtomicSpecies('S', 32.066, "S.pz-n-rrkjus_psl.0.1.UPF"),
[0.500000000, 0.288675130, 1.974192764],
)
AtomicPosition("S", [0.5, 0.28867513, 1.974192764], Bool[1, 1, 1])QuantumESPRESSOBase.PWscf.AtomicPositionsCard — TypeAtomicPositionsCard <: CardRepresent the ATOMIC_POSITIONS card in pw.x input files.
Arguments
data::AbstractVector{AtomicPosition}: A vector containingAtomicPositions.option::String="alat": allowed values are: "alat", "bohr", "angstrom", "crystal", and "crystal_sg".
QuantumESPRESSOBase.PWscf.CellParametersCard — TypeCellParametersCard(data::AbstractMatrix, option::String)Represent the CELL_PARAMETERS cards in PWscf and CP packages.
QuantumESPRESSOBase.PWscf.KMeshCard — TypeKMeshCard(data::MonkhorstPackGrid)Represent the K_POINTS card in Quantum ESPRESSO with Monkhorst–Pack grid.
QuantumESPRESSOBase.PWscf.GammaPointCard — TypeGammaPointCard()Represent the K_POINTS card in Quantum ESPRESSO with only Γ-point.
QuantumESPRESSOBase.PWscf.SpecialPointsCard — TypeSpecialPointsCard(data, option)Represent the K_POINTS card in Quantum ESPRESSO with a list of k-points.
Arguments
data::Vector{ReciprocalPoint}: a vector containingReciprocalPoints.option::String="tpiba": allowed values are: "tpiba", "automatic", "crystal", "gamma", "tpibab", "crystalb", "tpibac" and "crystalc".
QuantumESPRESSOBase.PWscf.PWInput — TypePWInput(control, system, electrons, ions, cell, atomic_species, atomic_positions, k_points, cell_parameters)Construct a PWInput which represents the input of program pw.x.
Arguments
control::ControlNamelist=ControlNamelist(): theCONTROLnamelist of the input. Optional.system::SystemNamelist=SystemNamelist(): theSYSTEMnamelist of the input. Optional.electrons::ElectronsNamelist=ElectronsNamelist(): theELECTRONSnamelist of the input. Optional.ions::IonsNamelist=IonsNamelist(): theIONSnamelist of the input. Optional.cell::CellNamelist=CellNamelist(): theCELLnamelist of the input. Optional.atomic_species::AtomicSpeciesCard: theATOMIC_SPECIEScard of the input. Must be provided explicitly.atomic_positions::AtomicPositionsCard: theATOMIC_POSITIONScard of the input. Must be provided explicitly.k_points::AbstractKPointsCard: theK_POINTScard of the input. Must be provided explicitly.cell_parameters::Union{Nothing,CellParametersCard}: theCELL_PARAMETERScard of the input. Must be eithernothingor aCellParametersCard.
Methods
QuantumESPRESSOBase.PWscf.isrequired — Functionisrequired(nml::Namelist)
isrequired(card::Card)Test whether a Namelist or a Card is required in PWInput.
QuantumESPRESSOBase.PWscf.isoptional — Functionisoptional(nml::Namelist)
isoptional(card::Card)Test whether a Namelist or a Card is optional in PWInput.
CrystallographyBase.cellvolume — Functioncellvolume(card)Return the cell volume of a CellParametersCard or RefCellParametersCard, in atomic unit.
It will throw an error if the option is "alat".
cellvolume(nml::SystemNamelist)Return the volume of the cell based on the information given in a SystemNamelist, in atomic unit.
cellvolume(input::PWInput)Return the volume of the cell based on the information given in a PWInput, in atomic unit.
QuantumESPRESSOBase.PWscf.convertoption — Functionconvertoption(card::AbstractCellParametersCard, new_option::AbstractString)Convert the option of an AbstractCellParametersCard from "bohr" to "angstrom", etc.
It does not support conversions between "alat" and others.
AbInitioSoftwareBase.getpseudodir — Functiongetpseudodir(nml::ControlNamelist)Get the directory storing the pseudopotential files.
getpseudodir(input::PWInput)Get the directory storing the pseudopotential files.
AbInitioSoftwareBase.listpotentials — Functionlistpotentials(card::AtomicSpeciesCard)List the pseudopotentials in an AtomicSpeciesCard.
listpotentials(input::PWInput)List the pseudopotentials in a PWInput.
QuantumESPRESSOBase.PWscf.getxmldir — Functiongetxmldir(nml::ControlNamelist)Return the path to directory storing the XML files.
QuantumESPRESSOBase.PWscf.listwfcfiles — Functionlistwfcfiles(nml::ControlNamelist, n=1)List all wave function files.