The Driver object

The Driver object is a simple way to call the QE/QEpy from python side. Here is how to define a simple scf job:

from qepy import Driver
driver = Driver('qe_in.in')

Here, the first argument specifies the name of QE input file.

List of all Methods

class qepy.driver.Driver(inputfile=None, comm=None, ldescf=False, iterative=False, task='scf', embed=None, prefix=None, outdir=None, logfile=None, qe_options=None, prog='pw', progress=False, atoms=None, needwf=True, **kwargs)[source]

The driver of QEpy.

Parameters
  • inputfile (str) -- Name of QE input file

  • comm (object) -- Parallel communicator

  • ldescf (bool) -- If True, also print the SCF correction term for each SCF cycle

  • iterative (bool) -- Run SCF or real-time TDDFT one iteration at a time

  • task (str) --

    Task to be performed by the driver :

    • 'scf' : Self consistent field

    • 'nscf' : initialization of Driver from previous SCF calculation

    • 'optical' : Optical absorption spectrum (real-time TDDFT with ce-tddft of Davide Ceresoli)

  • embed (object (Fortran)) -- embed object itialized in the Fortran side of QEpy needed to communicate with QE.

  • prefix (str) -- prefix of QE input/output

  • outdir (str) -- The output directory of QE

  • logfile (str, file or bool) --

    The screen output of QE. If it is a file descriptor, it should open with 'w+'.

  • outdir -- The output directory of QE

  • qe_options (dict) -- A dictionary with input parameters for QE to generate QE input file. See class QEInput

  • prog (str) -- The name of QE program, default is pw which is pw.x in QE.

  • progress (bool) -- If True, will continue run the QE without cleaning the workspace. Most times used for running TDDFT after scf.

  • atoms (object) -- An ase Atoms to generate the input file.

  • needwf (bool) -- If False, will not read wavefunctions and skip wavefunction-related initialization.

  • kwargs (dict) -- Other options

  • note:: (..) --

    • The bool gather parameter in functions means the data is gathered on root or distributed on all cpus.

    • Units in QEpy are Bohr and Rydberg. These differ from ASE's units. Please, be careful!

      • positions : Bohr

      • lattice : Bohr

      • energy : Ry

      • forces : Ry/Bohr

      • stress : Ry/Bohr**3

      • potential : Ry

      • density : 1.0/Bohr**3

calc_energy(**kwargs)[source]

Calculate the energy with PW2CASINO of QE. Use this only if you have a good reason to use it! For example embedding calculations.

check_convergence(**kwargs)[source]

Check the convergence of the SCF

create_array(gather=True, kind='rho')[source]

Returns an empty array in real space. Nota bene: this is for real-space arrays like the density (rho) and potentials.

classmethod data2field(data, cell=None, grid=None, rank=None)[source]

QE data to DFTpy DirectField. If data is None or small temporary array will return np.zeros(1).

diagonalize(print_level=2, nscf=False, **kwargs)[source]

Diagonalize the Hamiltonian

Parameters

print_level -- The level of output of QE

driver_initialize(**kwargs)[source]

Initialize the driver

Parameters
  • inputfile (str) -- Name of QE input file

  • commf (object) -- mpi4py parallel communicator to be sent to Fortran

  • task (str) --

    Task to be performed by the driver :

    • 'scf' : Self consistent field

    • 'nscf' : initialization of Driver from previous SCF calculation

    • 'optical' : Optical absorption spectrum (real-time TDDFT with ce-tddft of Davide Ceresoli)

  • qe_options (dict) -- A dictionary with input parameters for QE to generate QE input file.

electrons(original=False, **kwargs)[source]

Execute original QE routine "electrons"

end_scf(nscf=False, **kwargs)[source]

Ends the SCF and cleans the SCF workspace. Only run when in iterative mode (i.e., one cycle at a time)

classmethod field2data(field, data=None)[source]

DFTpy DirectField to QE data. If the input field is not 3d array, will return np.zeros((1, 1)).

classmethod get_ase_atoms()[source]

Return the atom.Atoms from QE to ASE (using ASE units).

classmethod get_bz_k_points()[source]

Return all the k-points in the 1. Brillouin zone.

The coordinates are relative to reciprocal latice vectors.

get_core_density(gather=True, out=None)[source]

Returns core density array in real space.

get_density(gather=True, out=None)[source]

Returns valence pseudodensity array in real space.

get_density_functional(gather=True, out=None, add=False, **kwargs)[source]

Return effective potential information.

Note :

Then final potential is saved in the v_obj

classmethod get_dftpy_grid(nr=None, cell=None, mp=None, **kwargs)[source]

Return the DFTpy DirectGrid from QE.

classmethod get_dftpy_ions()[source]

Return the DFTpy.Ions from QE. (with DFTpy units)

get_dipole_tddft()[source]

Return the total dipole computed by a TDDFT run. Nota bene: only available during TDDFT runs.

classmethod get_ecutrho()[source]

Return the cutoff for density.

classmethod get_eigenvalues(kpt=0, spin=0)[source]

Return eigenvalue array.

get_elf(gather=True, out=None, **kwargs)[source]

Return electron localization function.

get_energy(**kwargs)[source]

Return the total energy. Nota bene: Only use for regular QE runs (i.e., not when doing embedding).

get_exchange_correlation(gather=True, out=None, tau=None, **kwargs)[source]

Return information about the exchange-correlation: LDA, GGA: tuple (potential, energy, v*rho) MGGA: tuple (potential, energy, v*rho, tau)

TODO :

The interface will changed in the new version of QE!

Note :

For metaGGA, only return scattered tau

classmethod get_fermi_level()[source]

Return the Fermi level.

get_forces(icalc=0, ignore=(), **kwargs)[source]

Return the total forces. (n, 3)

Parameters
  • icalc (int) --

    icalc

    without

    bin

    0

    all

    000

    1

    no ewald

    001

    2

    no localpp

    010

    4

    no nlcc

    100

  • ignore (list) --

    ignore some forces, which does same job as icalc.

    • ewald (1)

    • localpp (2)

    • nlcc (4)

get_hartree(gather=True, out=None, add=False, **kwargs)[source]

Return information about Hartree energy: tuple (potential, energy, total charge)

classmethod get_ibz_k_points()[source]

Return k-points in the irreducible part of the Brillouin zone.

The coordinates are relative to reciprocal latice vectors.

classmethod get_ions_lattice()[source]

Return the matrix of the lattice vectors in Bohrs.

classmethod get_ions_positions()[source]

Return the cartesian positions of ions in Bohrs.

classmethod get_ions_symbols()[source]

Return the symbols of ions.

classmethod get_k_point_weights()[source]

Weights of the k-points.

The sum of all weights is one.

get_kinetic_energy_density(gather=True, out=None)[source]

Returns KS kinetic energy density array in real space.

get_local_pp(gather=True, out=None, **kwargs)[source]

Return local component of the pseudopotential.

classmethod get_magnetic_moment(**kwargs)[source]

Return the total magnetic moment.

classmethod get_number_of_bands()[source]

Return the number of bands.

classmethod get_number_of_grid_points(gather=True)[source]

Return the shape of arrays.

classmethod get_number_of_k_points()[source]

Return the number of kpoints.

classmethod get_number_of_spins()[source]

Return the number of spins in the calculation.

Spin-paired calculations: 1, spin-polarized calculation: 2.

classmethod get_occupation_numbers(kpt=0, spin=0)[source]

Return occupation number array.

get_output()[source]

Return the output of QE.

It depends on the logfile of the initialization of the driver :

  • None : return None.

  • str : return all outputs.

  • True : It will return the output from last time.

get_potential_energy(**kwargs)[source]

Return the potential energy.

get_pseudo_density(spin=None, pad=True, gather=True)[source]

Return pseudo-density array.

If spin is not given, then the total density is returned. Otherwise, the spin up or down density is returned (spin=0 or 1).

classmethod get_pseudo_wave_function(band=None, kpt=0, spin=0, broadcast=True, pad=True)[source]

Return pseudo-wave-function array.

get_rdg(gather=True, out=None, **kwargs)[source]

Return reduced density gradient.

get_scf_error(**kwargs)[source]

Return the error of the SCF compared to previous cycle

get_scf_steps(**kwargs)[source]

Return the number of SCF steps

classmethod get_spin_polarized()[source]

Is it a spin-polarized calculation?

classmethod get_stress(**kwargs)[source]

Return the stress (3, 3).

classmethod get_volume()[source]

Return the volume.

get_wave_function(band=None, kpt=0)[source]

Returns wave-function array in real space.

classmethod get_xc_functional()[source]

Return the XC-functional identifier.

'LDA', 'PBE', ...

property is_root

Whether it is the root.

If the comm is set, root is rank == 0. If comm not set, root is ionode of QE.

mix(mix_coef=0.7, print_level=2)[source]

Mix the density with the QE density mixing

Parameters

print_level -- The level of output of QE

non_scf(**kwargs)[source]

Single "non-selfconsistent" calculation

pwscf_restart(oldxml=False, starting_pot='file', starting_wfc='file')[source]

Read PW ouput/restart files.

Parameters

oldxml (bool) --

  • True : read the old format xml file (qe<6.4)

  • False : read the new format xml file (qe>6.3)

restart(prog=None, **kwargs)[source]

Restart the driver losing all information about the previous driver

save(what='all', **kwargs)[source]

Save the QE data to the disk in original QE files that can be accessed via post-processing or read with a QEpy driver

Parameters

what (str) --

  • what = 'all'write xml data file, charge density, wavefunctions

    (for final data);

  • what = 'config'write xml data file and charge density; also,

    for nks=1, wavefunctions in plain binary format (see why in comments below). For intermediate or incomplete results;

  • what = 'config-nowf'write xml data file and charge density only;

    (save density);

  • what = 'config-init'write xml data file only excluding final results

    (for dry run, can be called at early stages).

see PW/src/punch.f90

scf(print_level=2, maxiter=None, original=False, nscf=False, **kwargs)[source]

Run the scf/tddft until converged or reached the maximum number of iterations

classmethod set_density(density, gather=True, **kwargs)[source]

Set density array in real space.

set_external_potential(potential, exttype=None, gather=True, extene=None, **kwargs)[source]

Set an external potential in addition to the ones already included in the QEpy run according to the logic enumerated below.

Parameters
  • potential ((nnr, nspin)) -- The external potential

  • exttype (list or int) --

    The type of external potential. It can be a list of name or a integer. e.g. exttype = ('localpp', 'xc') or exttype = 5. external stands for an additional external potential.

    type

    potential

    bin

    0

    external

    000

    1

    localpp

    001

    2

    hartree

    010

    4

    xc

    100

stop(exit_status=0, what='all', print_flag=0, **kwargs)[source]

Stop the driver. This must be done anytime a new driver is created. This method is invoked automatically if a running driver is detected. Only one driver can run at any given time.

Parameters
  • exit_status (int) -- 0 : QE will remove it temporary files

  • print_flag (int) -- 0 : Not print time informations

  • what (str) -- see qepy.driver.Driver.save().

classmethod sum_band(occupations=None, **kwargs)[source]

Same as sum_band of QE with input occupations:

Parameters

occupations (np.ndarray (nbnd, nk)) -- occupation numbers

classmethod switch_nlpp(nhm=0, nbetam=0, nkb=0, nh=None, **kwargs)[source]

Switch on/off the nonlocal part of the pseudopotential

tddft_initialize(inputfile=None, commf=None, **kwargs)[source]

Initialize the tddft

Parameters
  • inputfile (str) -- Name of QE input file, which also contains &inputtddft section.

  • commf (object) -- mpi4py parallel communicator to be sent to Fortran

tddft_restart(istep=None, **kwargs)[source]

Restarts the TDDFT from previous interrupted run.

Parameters

istep (int) -- Start number of steps, just for output.

tddft_stop(exit_status=0, what='no', print_flag=0, **kwargs)[source]

Stops TDDFT run

classmethod update_exchange_correlation(xc=None, libxc=None, **kwargs)[source]

Switch XC on the fly

update_ions(positions=None, lattice=None, update=0, **kwargs)[source]

update the ions of QE. This is used for dynamic runs.

Parameters
  • positions (np.ndarray (n, 3)) -- positions of ions

  • lattice (3x3 matrix) -- lattice of ions

  • update (int) --

    • 0 : update everything (default)

    • 1 : only update atomic information, not update potential