FHI-aims
Introduction
FHI-aims is a all-electron full-potential density functional theory code using a numeric local orbital basis set.
Running the Calculator
The default initialization command for the FHI-aims calculator is
- class ase.calculators.aims.Aims(profile=None, directory='.', **kwargs)[source]
 Construct the FHI-aims calculator.
The keyword arguments (kwargs) can be one of the ASE standard keywords: ‘xc’, ‘kpts’ and ‘smearing’ or any of FHI-aims’ native keywords.
Arguments:
- cubes: AimsCube object
 Cube file specification.
- tier: int or array of ints
 Set basis set tier for all atomic species.
- plus_udict
 For DFT+U. Adds a +U term to one specific shell of the species.
- kwargsdict
 Any of the base class arguments.
In order to run a calculation, you have to ensure that at least the
following str variables are specified, either in the initialization
or as shell environment variables:
keyword  | 
description  | 
|---|---|
  | 
The full command required to run FHI-aims from
a shell, including anything to do with an MPI
wrapper script and the number of tasks, e.g.:
  | 
  | 
Directory where the species are located, e.g.:
  | 
  | 
which exchange-correlation functional is used.  | 
List of keywords
This is a non-exclusive list of keywords for the control.in file
that can be addresses from within ASE. The meaning for these keywords is
exactly the same as in FHI-aims, please refer to its manual for help on
their use.
One thing that should be mentioned is that keywords with more than
one option have been implemented as tuples/lists, eg.
k_grid=(12,12,12) or relativistic=('atomic_zora','scalar').
In those cases, specifying a single string containing all the options is also possible.
None of the keywords have any default within ASE, but do check the defaults set by FHI-aims.
Example keywords describing the computational method used:
keyword  | 
type  | 
|---|---|
  | 
str  | 
  | 
float  | 
  | 
str  | 
  | 
list  | 
  | 
bool  | 
  | 
str  | 
  | 
list  | 
Note
Any argument can be changed after the initial construction of the calculator, simply by setting it with the method
>>> calc.set(keyword=value)
Volumetric Data Output
The class
- class ase.calculators.aims.AimsCube(origin=(0, 0, 0), edges=[(0.1, 0.0, 0.0), (0.0, 0.1, 0.0), (0.0, 0.0, 0.1)], points=(50, 50, 50), plots=())[source]
 Object to ensure the output of cube files, can be attached to Aims object
parameters:
- origin, edges, points:
 Same as in the FHI-aims output
- plots:
 what to print, same names as in FHI-aims
describes an object that takes care of the volumetric output requests within FHI-aims. An object of this type can be attached to the main Aims() object as an option.
The possible arguments for AimsCube are:
keyword  | 
type  | 
|---|---|
  | 
list  | 
  | 
3x3-array  | 
  | 
list  | 
  | 
list  | 
The possible values for the entry of plots are discussed in detail in the FHI-aims manual, see below for an example.
Example
Here is an example of how to obtain the geometry of a water molecule,
assuming ASE_AIMS_COMMAND and AIMS_SPECIES_DIR are set:
:git:`ase/test/calculator/aims/test_H2O_aims.py`.
import pytest
from ase import Atoms
from ase.calculators.aims import AimsCube
from ase.optimize import QuasiNewton
@pytest.mark.calculator('aims')
def test_H2O_aims(factory):
    water = Atoms('HOH', [(1, 0, 0), (0, 0, 0), (0, 1, 0)])
    water_cube = AimsCube(points=(29, 29, 29),
                          plots=('total_density',
                                 'delta_density',
                                 'eigenstate 5',
                                 'eigenstate 6'))
    calc = factory.calc(
        xc='LDA',
        output=['dipole'],
        sc_accuracy_etot=1e-2,
        sc_accuracy_eev=1e-1,
        sc_accuracy_rho=1e-2,
        sc_accuracy_forces=1e-1,
        cubes=water_cube
    )
    water.calc = calc
    dynamics = QuasiNewton(water)
    dynamics.run(fmax=0.2)