QuantumESPRESSOBase.PWscf module

Types

CrystallographyBase.LatticeType
Lattice(nml::SystemNamelist)

Create a Lattice from a SystemNamelist.

source
Lattice(card::CellParametersCard)

Create a Lattice from a CellParametersCard.

source
Lattice(card::PWInput)

Create a Lattice from a PWInput.

source
QuantumESPRESSOBase.PWscf.ControlNamelistType
ControlNamelist(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.

source
QuantumESPRESSOBase.PWscf.SystemNamelistType
SystemNamelist(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.

source
QuantumESPRESSOBase.PWscf.ElectronsNamelistType
ElectronsNamelist(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.

source
QuantumESPRESSOBase.PWscf.IonsNamelistType
IonsNamelist(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".

source
QuantumESPRESSOBase.PWscf.CellNamelistType
CellNamelist(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".

source
QuantumESPRESSOBase.PWscf.DosNamelistType
DosNamelist(prefix, outdir, ngauss, degauss, Emin, Emax, DeltaE, fildos)
DosNamelist(; kwargs...)
DosNamelist(::DosNamelist; kwargs...)

Represent the DOS namelist of dos.x.

source
QuantumESPRESSOBase.PWscf.BandsNamelistType
BandsNamelist(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.

source
QuantumESPRESSOBase.PWscf.AtomicSpeciesType
AtomicSpecies(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")
source
QuantumESPRESSOBase.PWscf.AtomicPositionType
AtomicPosition(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])
source
QuantumESPRESSOBase.PWscf.AtomicPositionsCardType
AtomicPositionsCard <: Card

Represent the ATOMIC_POSITIONS card in pw.x input files.

Arguments

  • data::AbstractVector{AtomicPosition}: A vector containing AtomicPositions.
  • option::String="alat": allowed values are: "alat", "bohr", "angstrom", "crystal", and "crystal_sg".
source
QuantumESPRESSOBase.PWscf.SpecialPointsCardType
SpecialPointsCard(data, option)

Represent the K_POINTS card in Quantum ESPRESSO with a list of k-points.

Arguments

  • data::Vector{ReciprocalPoint}: a vector containing ReciprocalPoints.
  • option::String="tpiba": allowed values are: "tpiba", "automatic", "crystal", "gamma", "tpibab", "crystalb", "tpibac" and "crystalc".
source
QuantumESPRESSOBase.PWscf.PWInputType
PWInput(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(): the CONTROL namelist of the input. Optional.
  • system::SystemNamelist=SystemNamelist(): the SYSTEM namelist of the input. Optional.
  • electrons::ElectronsNamelist=ElectronsNamelist(): the ELECTRONS namelist of the input. Optional.
  • ions::IonsNamelist=IonsNamelist(): the IONS namelist of the input. Optional.
  • cell::CellNamelist=CellNamelist(): the CELL namelist of the input. Optional.
  • atomic_species::AtomicSpeciesCard: the ATOMIC_SPECIES card of the input. Must be provided explicitly.
  • atomic_positions::AtomicPositionsCard: the ATOMIC_POSITIONS card of the input. Must be provided explicitly.
  • k_points::AbstractKPointsCard: the K_POINTS card of the input. Must be provided explicitly.
  • cell_parameters::Union{Nothing,CellParametersCard}: the CELL_PARAMETERS card of the input. Must be either nothing or a CellParametersCard.
source

Methods

CrystallographyBase.cellvolumeFunction
cellvolume(card)

Return the cell volume of a CellParametersCard or RefCellParametersCard, in atomic unit.

Warning

It will throw an error if the option is "alat".

source
cellvolume(nml::SystemNamelist)

Return the volume of the cell based on the information given in a SystemNamelist, in atomic unit.

source
cellvolume(input::PWInput)

Return the volume of the cell based on the information given in a PWInput, in atomic unit.

source
QuantumESPRESSOBase.PWscf.convertoptionFunction
convertoption(card::AbstractCellParametersCard, new_option::AbstractString)

Convert the option of an AbstractCellParametersCard from "bohr" to "angstrom", etc.

Warning

It does not support conversions between "alat" and others.

source