The GPAW calculator object

class gpaw.calculator.GPAW(restart=None, *, label=None, timer=None, communicator=None, txt='?', parallel=None, **kwargs)[source]

This is the ASE-calculator frontend for doing a GPAW calculation.

Basic calculator implementation.

restart: str

Prefix for restart file. May contain a directory. Default is None: don’t restart.

ignore_bad_restart_file: bool

Deprecated, please do not use. Passing more than one positional argument to Calculator() is deprecated and will stop working in the future. Ignore broken or missing restart file. By default, it is an error if the restart file is missing or broken.

directory: str or PurePath

Working directory in which to read and write files and perform calculations.

label: str

Name used for all files. Not supported by all calculators. May contain a directory, but please use the directory parameter for that instead.

atoms: Atoms object

Optional Atoms object to which the calculator will be attached. When restarting, atoms will get its positions and unit-cell updated from file.

add_wannier_correction(Z_nn, G_c, u, u1, nbands=None)[source]

Calculate the correction to the wannier integrals.

See: (Eq. 27 ref1):

              -i G.r
Z   = <psi | e      |psi >
 nm       n             m

               __                __
       ~      \              a  \     a*   a    a
Z    = Z    +  ) exp[-i G . R ]  )   P   dO    P
 nmx    nmx   /__            x  /__   ni   ii'  mi'

               a                 ii'

Note that this correction is an approximation that assumes the exponential varies slowly over the extent of the augmentation sphere.

ref1: Thygesen et al, Phys. Rev. B 72, 125119 (2005)

attach(function, n=1, *args, **kwargs)[source]

Register observer function to run during the SCF cycle.

Call function using args and kwargs as arguments.

If n is positive, then function will be called every n SCF iterations + the final iteration if it would not be otherwise

If n is negative, then function will only be called on iteration abs(n).

If n is 0, then function will only be called on convergence

band_structure()

Create band-structure object for plotting.

calculate(atoms=None, properties=['energy'], system_changes=['cell'])[source]

Do the calculation.

properties: list of str

List of what needs to be calculated. Can be any combination of ‘energy’, ‘forces’, ‘stress’, ‘dipole’, ‘charges’, ‘magmom’ and ‘magmoms’.

system_changes: list of str

List of what has changed since last calculation. Can be any combination of these six: ‘positions’, ‘numbers’, ‘cell’, ‘pbc’, ‘initial_charges’ and ‘initial_magmoms’.

Subclasses need to implement this, but can ignore properties and system_changes if they want. Calculated properties should be inserted into results dictionary like shown in this dummy example:

self.results = {'energy': 0.0,
                'forces': np.zeros((len(atoms), 3)),
                'stress': np.zeros(6),
                'dipole': np.zeros(3),
                'charges': np.zeros(len(atoms)),
                'magmom': 0.0,
                'magmoms': np.zeros(len(atoms))}

The subclass implementation should first call this implementation to set the atoms attribute and create any missing directories.

calculate_numerical_forces(atoms, d=0.001)

Calculate numerical forces using finite difference.

All atoms will be displaced by +d and -d in all directions.

calculate_numerical_stress(atoms, d=1e-06, voigt=True)

Calculate numerical stress using finite difference.

calculate_properties(atoms, properties)

This method is experimental; currently for internal use.

call_observers(iter, final=False)[source]

Call all registered callback functions.

check_state(atoms, tol=1e-12)[source]

Check for any system changes since last calculation.

converge_wave_functions()[source]

Converge the wave-functions if not present.

default_parameters: Dict[str, Any] = {'background_charge': None, 'basis': {}, 'charge': 0, 'convergence': {'bands': 'occupied', 'density': 0.0001, 'eigenstates': 4e-08, 'energy': 0.0005}, 'dtype': None, 'eigensolver': None, 'experimental': {'kpt_refine': None, 'magmoms': None, 'niter_fixdensity': 0, 'reuse_wfs_method': 'paw', 'soc': None}, 'external': None, 'filter': None, 'fixdensity': False, 'gpts': None, 'h': None, 'hund': False, 'kpts': [(0.0, 0.0, 0.0)], 'maxiter': 333, 'mixer': None, 'mode': None, 'nbands': None, 'occupations': None, 'poissonsolver': None, 'random': False, 'setups': {}, 'spinpol': None, 'symmetry': {'do_not_symmetrize_the_density': None, 'point_group': True, 'symmorphic': True, 'time_reversal': True, 'tolerance': 1e-07}, 'verbose': 0, 'xc': 'LDA'}

Default parameters

discard_results_on_any_change = False

Whether we purge the results following any change in the set() method.

dos(soc=False, theta=0.0, phi=0.0, shift_fermi_level=True)[source]

Create DOS-calculator.

Default is to shift_fermi_level to 0.0 eV. For soc=True, angles can be given in degrees.

embed(q_p, rc=0.2, rc2=inf, width=1.0)[source]

Embed QM region in point-charges.

estimate_memory(mem)[source]

Estimate memory use of this object.

fixed_density(*, update_fermi_level=False, communicator=None, txt='-', parallel=None, **kwargs)[source]

Create new calculator and do SCF calculation with fixed density.

Returns a new GPAW object fully converged.

Useful for band-structure calculations. Given a ground-state calculation, gs_calc, one can do:

bs_calc = gs_calc.fixed_density(kpts=<path>,
                                symmetry='off')
bs = bs_calc.get_band_structure()
get_all_electron_density(spin=None, gridrefinement=2, pad=True, broadcast=True, collect=True, skip_core=False)[source]

Return reconstructed all-electron density array.

get_all_electron_ldos(mol, spin=0, npts=201, width=None, wf_k=None, P_aui=None, lc=None, raw=False)[source]

The Projected Density of States, using all-electron wavefunctions.

Projects onto a pseudo_wavefunctions (wf_k) corresponding to some band n and uses P_aui ([paw.nuclei[a].P_uni[:,n,:] for a in atoms]) to obtain the all-electron overlaps. Instead of projecting onto a wavefunction, a molecular orbital can be specified by a linear combination of weights (lc)

get_atomic_electrostatic_potentials()[source]

Return the electrostatic potential at the atomic sites.

Return list of energies in eV, one for each atom:

\[Y_{00} \int d\mathbf{r} \tilde{v}_{H} (\mathbf{r}) \hat{g}^{a}_{00} (\mathbf{r}-\mathbf{R}^{a} )\]
get_bz_k_points()[source]

Return the k-points.

get_bz_to_ibz_map()[source]

Return indices from BZ to IBZ.

get_dos(spin=0, npts=201, width=None)[source]

The total DOS.

Fold eigenvalues with Gaussians, and put on an energy grid.

returns an (energies, dos) tuple, where energies are relative to the vacuum level for non-periodic systems, and the average potential for periodic systems.

get_effective_potential(spin=0, pad=True, broadcast=True)[source]

Return pseudo effective-potential.

get_eigenvalues(kpt=0, spin=0, broadcast=True)[source]

Return eigenvalue array.

get_electrostatic_corrections()[source]

Calculate PAW correction to average electrostatic potential.

get_electrostatic_potential()[source]

Return the electrostatic potential.

This is the potential from the pseudo electron density and the PAW-compensation charges. So, the electrostatic potential will only be correct outside the PAW augmentation spheres.

get_fermi_level()[source]

Return the Fermi-level.

get_fermi_levels()[source]

Return the Fermi-levels in case of fixed-magmom.

get_homo_lumo(spin=None)[source]

Return HOMO and LUMO eigenvalues.

By default, return the true HOMO-LUMO eigenvalues (spin=None).

If spin is 0 or 1, return HOMO-LUMO eigenvalues taken among only those states with the given spin.

get_ibz_k_points()[source]

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

get_k_point_weights()[source]

Weights of the k-points.

The sum of all weights is one.

get_lcao_dos(atom_indices=None, basis_indices=None, npts=201, width=None)[source]

Get density of states projected onto orbitals in LCAO mode.

basis_indices is a list of indices of basis functions on which to project. To specify all basis functions on a set of atoms, you can supply atom_indices instead. Both cannot be given simultaneously.

get_magnetic_moments(atoms=None)

Calculate magnetic moments projected onto atoms.

get_number_of_bands()[source]

Return the number of bands.

get_occupation_numbers(kpt=0, spin=0, broadcast=True)[source]

Return occupation array.

get_orbital_ldos(a, spin=0, angular='spdf', npts=201, width=None, nbands=None, spinorbit=False)[source]

The Local Density of States, using atomic orbital basis functions.

Project wave functions onto an atom orbital at atom a, and use this as weight when summing the eigenvalues.

The atomic orbital has angular momentum angular, which can be ‘s’, ‘p’, ‘d’, ‘f’, or any combination (e.g. ‘sdf’).

An integer value for angular can also be used to specify a specific projector function to project onto.

Setting nbands limits the number of bands included. This speeds up the calculation if one has many bands in the calculator but is only interested in the DOS at low energies.

get_projections(locfun, spin=0)[source]

Project wave functions onto localized functions

Determine the projections of the Kohn-Sham eigenstates onto specified localized functions of the format:

locfun = [[spos_c, l, sigma], [...]]

spos_c can be an atom index, or a scaled position vector. l is the angular momentum, and sigma is the (half-) width of the radial gaussian.

Return format is:

f_kni = <psi_kn | f_i>

where psi_kn are the wave functions, and f_i are the specified localized functions.

As a special case, locfun can be the string ‘projectors’, in which case the bound state projectors are used as localized functions.

get_property(name, atoms=None, allow_calculation=True)

Get the named property.

get_pseudo_density(spin=None, gridrefinement=1, pad=True, broadcast=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).

get_pseudo_density_corrections()[source]

Integrated density corrections.

Returns the integrated value of the difference between the pseudo- and the all-electron densities at each atom. These are the numbers you should add to the result of doing e.g. Bader analysis on the pseudo density.

get_pseudo_valence_density(spin=None, gridrefinement=1, pad=True, broadcast=True)

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).

get_pseudo_wave_function(band=0, kpt=0, spin=0, broadcast=True, pad=True, periodic=False)[source]

Return pseudo-wave-function array.

Units: 1/Angstrom^(3/2)

get_spin_polarized()[source]

Is it a spin-polarized calculation?

get_stresses(atoms=None)

the calculator should return intensive stresses, i.e., such that stresses.sum(axis=0) == stress

get_wannier_integrals(s, k, k1, G_c, nbands=None)[source]

Calculate integrals for maximally localized Wannier functions.

get_wannier_localization_matrix(nbands, dirG, kpoint, nextkpoint, G_I, spin)[source]

Calculate integrals for maximally localized Wannier functions.

get_wigner_seitz_densities(spin)[source]

Get the weight of the spin-density in Wigner-Seitz cells around each atom.

The density assigned to each atom is relative to the neutral atom, i.e. the density sums to zero.

get_wigner_seitz_ldos(a, spin=0, npts=201, width=None)[source]

The Local Density of States, using a Wigner-Seitz basis function.

Project wave functions onto a Wigner-Seitz box at atom a, and use this as weight when summing the eigenvalues.

get_xc_functional()[source]

Return the XC-functional identifier.

‘LDA’, ‘PBE’, …

icalculate(atoms=None, properties=['energy'], system_changes=['cell'])[source]

Calculate things.

ignored_changes: Set[str] = {}

Properties of Atoms which we ignore for the purposes of cache

implemented_properties: List[str] = ['energy', 'free_energy', 'forces', 'stress', 'dipole', 'magmom', 'magmoms']

Properties calculator can handle (energy, forces, …)

initial_wannier(initialwannier, kpointgrid, fixedstates, edf, spin, nbands)[source]

Initial guess for the shape of wannier functions.

Use initial guess for wannier orbitals to determine rotation matrices U and C.

initialize(atoms=None, reading=False)[source]

Inexpensive initialization.

initialize_positions(atoms=None)[source]

Update the positions of the atoms.

linearize_to_xc(newxc)[source]

Linearize Hamiltonian to difference XC functional.

Used in real time TDDFT to perform calculations with various kernels.

new(timer=None, communicator=None, txt='-', parallel=None, **kwargs)[source]

Create a new calculator, inheriting input parameters.

The txt file and timer are the only input parameters to be created anew. Internal variables, such as the density or the wave functions, are not reused either.

For example, to perform an identical calculation with a parameter changed (e.g. changing XC functional to PBE):

new_calc = calc.new(xc='PBE')
atoms.calc = new_calc
print_memory_estimate(log=None, maxdepth=-1)[source]

Print estimated memory usage for PAW object and components.

maxdepth is the maximum nesting level of displayed components.

The PAW object must be initialize()’d, but needs not have large arrays allocated.

read(filename)[source]

Read atoms, parameters and calculated properties from output file.

Read result from self.label file. Raise ReadError if the file is not there. If the file is corrupted or contains an error message from the calculation, a ReadError should also be raised. In case of succes, these attributes must set:

atoms: Atoms object

The state of the atoms from last calculation.

parameters: Parameters object

The parameter dictionary.

results: dict

Calculated properties like energy and forces.

The FileIOCalculator.read() method will typically read atoms and parameters and get the results dict by calling the read_results() method.

reset()

Clear all information from old calculation.

set(**kwargs)[source]

Change parameters for calculator.

Examples:

calc.set(xc='PBE')
calc.set(nbands=20, kpts=(4, 1, 1))
set_label(label)

Set label and convert label to directory and prefix.

Examples:

  • label=’abc’: (directory=’.’, prefix=’abc’)

  • label=’dir1/abc’: (directory=’dir1’, prefix=’abc’)

  • label=None: (directory=’.’, prefix=None)

set_positions(atoms=None)[source]

Update the positions of the atoms and initialize wave functions.

write(filename, mode='')[source]

Write calculator object to a file.

Parameters:
  • filename – File to be written

  • mode – Write mode. Use mode='all' to include wave functions in the file.