Timings and parallelization

This section summarizes the options DFTK offers to monitor and influence performance of the code.

Timing measurements

By default DFTK uses TimerOutputs.jl to record timings, memory allocations and the number of calls for selected routines inside the code. These numbers are accessible in the object DFTK.timer. Since the timings are automatically accumulated inside this datastructure, any timing measurement should first reset this timer before running the calculation of interest.

For example to measure the timing of an SCF:

DFTK.reset_timer!(DFTK.timer)
scfres = self_consistent_field(basis, tol=1e-5)

DFTK.timer
────────────────────────────────────────────────────────────────────────────────
                                       Time                    Allocations      
                              ───────────────────────   ────────────────────────
      Tot / % measured:            287ms /  50.7%           68.4MiB /  80.6%    

Section               ncalls     time    %tot     avg     alloc    %tot      avg
────────────────────────────────────────────────────────────────────────────────
self_consistent_field      1    146ms  100.0%   146ms   55.1MiB  100.0%  55.1MiB
  LOBPCG                  27   67.3ms   46.2%  2.49ms   15.6MiB   28.2%   590KiB
    DftHamiltonian...     75   51.7ms   35.5%   689μs   4.49MiB    8.1%  61.3KiB
      local+kinetic      491   49.7ms   34.2%   101μs    230KiB    0.4%     480B
      nonlocal            75    878μs    0.6%  11.7μs    317KiB    0.6%  4.22KiB
    rayleigh_ritz         48   4.21ms    2.9%  87.7μs   1.72MiB    3.1%  36.8KiB
      ortho!              48    415μs    0.3%  8.64μs    252KiB    0.4%  5.24KiB
    ortho! X vs Y         69   3.97ms    2.7%  57.6μs   1.87MiB    3.4%  27.8KiB
      ortho!             135   1.74ms    1.2%  12.9μs   1.00MiB    1.8%  7.58KiB
    Update residuals      75    790μs    0.5%  10.5μs   1.14MiB    2.1%  15.5KiB
    preconditioning       75    784μs    0.5%  10.5μs    224KiB    0.4%  2.98KiB
    ortho!                27    435μs    0.3%  16.1μs    160KiB    0.3%  5.91KiB
  compute_density          9   47.2ms   32.4%  5.24ms   6.25MiB   11.3%   711KiB
    symmetrize_ρ           9   41.3ms   28.4%  4.59ms   5.01MiB    9.1%   569KiB
  energy_hamiltonian      10   14.2ms    9.8%  1.42ms   14.3MiB   26.0%  1.43MiB
    ene_ops               10   12.7ms    8.7%  1.27ms   10.5MiB   19.0%  1.05MiB
      ene_ops: xc         10   10.5ms    7.2%  1.05ms   4.68MiB    8.5%   479KiB
      ene_ops: har...     10   1.45ms    1.0%   145μs   4.51MiB    8.2%   462KiB
      ene_ops: non...     10    239μs    0.2%  23.9μs    148KiB    0.3%  14.8KiB
      ene_ops: kin...     10    189μs    0.1%  18.9μs   96.9KiB    0.2%  9.69KiB
      ene_ops: local      10    130μs    0.1%  13.0μs    938KiB    1.7%  93.8KiB
  energy                   9   11.8ms    8.1%  1.31ms   9.43MiB   17.1%  1.05MiB
    ene_ops: xc            9   9.65ms    6.6%  1.07ms   4.21MiB    7.6%   479KiB
    ene_ops: hartree       9   1.35ms    0.9%   150μs   4.06MiB    7.4%   462KiB
    ene_ops: nonlocal      9    292μs    0.2%  32.5μs    148KiB    0.3%  16.5KiB
    ene_ops: kinetic       9    202μs    0.1%  22.4μs   91.1KiB    0.2%  10.1KiB
    ene_ops: local         9    121μs    0.1%  13.4μs    845KiB    1.5%  93.8KiB
  Anderson acceler...      8   2.21ms    1.5%   277μs   4.45MiB    8.1%   570KiB
  ortho_qr                 3    139μs    0.1%  46.5μs    101KiB    0.2%  33.7KiB
  χ0Mixing                 9   41.6μs    0.0%  4.62μs   41.1KiB    0.1%  4.56KiB
enforce_real!              1   25.0μs    0.0%  25.0μs   1.69KiB    0.0%  1.69KiB
────────────────────────────────────────────────────────────────────────────────

The output produced when printing or displaying the DFTK.timer now shows a nice table summarising total time and allocations as well as a breakdown over individual routines.

Timing measurements and stack traces

Timing measurements have the unfortunate disadvantage that they alter the way stack traces look making it sometimes harder to find errors when debugging. For this reason timing measurements can be disabled completely (i.e. not even compiled into the code) by setting the package-level preference DFTK.set_timer_enabled!(false). You will need to restart your Julia session afterwards to take this into account.

Rough timing estimates

A very (very) rough estimate of the time per SCF step (in seconds) can be obtained with the following function. The function assumes that FFTs are the limiting operation and that no parallelisation is employed.

function estimate_time_per_scf_step(basis::PlaneWaveBasis)
    # Super rough figure from various tests on cluster, laptops, ... on a 128^3 FFT grid.
    time_per_FFT_per_grid_point = 30 #= ms =# / 1000 / 128^3

    (time_per_FFT_per_grid_point
     * prod(basis.fft_size)
     * length(basis.kpoints)
     * div(basis.model.n_electrons, DFTK.filled_occupation(basis.model), RoundUp)
     * 8  # mean number of FFT steps per state per k-point per iteration
     )
end

"Time per SCF (s):  $(estimate_time_per_scf_step(basis))"
"Time per SCF (s):  0.008009033203124998"

Options for parallelization

At the moment DFTK offers two ways to parallelize a calculation, firstly shared-memory parallelism using threading and secondly multiprocessing using MPI (via the MPI.jl Julia interface). MPI-based parallelism is currently only over $k$-points, such that it cannot be used for calculations with only a single $k$-point. Otherwise combining both forms of parallelism is possible as well.

The scaling of both forms of parallelism for a number of test cases is demonstrated in the following figure. These values were obtained using DFTK version 0.1.17 and Julia 1.6 and the precise scalings will likely be different depending on architecture, DFTK or Julia version. The rough trends should, however, be similar.

The MPI-based parallelization strategy clearly shows a superior scaling and should be preferred if available.

MPI-based parallelism

Currently DFTK uses MPI to distribute on $k$-points only. This implies that calculations with only a single $k$-point cannot use make use of this. For details on setting up and configuring MPI with Julia see the MPI.jl documentation.

  1. First disable all threading inside DFTK, by adding the following to your script running the DFTK calculation:

    using DFTK
    disable_threading()
  2. Run Julia in parallel using the mpiexecjl wrapper script from MPI.jl:

    mpiexecjl -np 16 julia myscript.jl

    In this -np 16 tells MPI to use 16 processes and -t 1 tells Julia to use one thread only. Notice that we use mpiexecjl to automatically select the mpiexec compatible with the MPI version used by MPI.jl.

As usual with MPI printing will be garbled. You can use

DFTK.mpi_master() || (redirect_stdout(); redirect_stderr())

at the top of your script to disable printing on all processes but one.

MPI-based parallelism not fully supported

While standard procedures (such as the SCF or band structure calculations) fully support MPI, not all routines of DFTK are compatible with MPI yet and will throw an error when being called in an MPI-parallel run. In most cases there is no intrinsic limitation it just has not yet been implemented. If you require MPI in one of our routines, where this is not yet supported, feel free to open an issue on github or otherwise get in touch.

Thread-based parallelism

Threading in DFTK currently happens on multiple layers distributing the workload over different $k$-points, bands or within an FFT or BLAS call between threads. At its current stage our scaling for thread-based parallelism is worse compared MPI-based and therefore the parallelism described here should only be used if no other option exists. To use thread-based parallelism proceed as follows:

  1. Ensure that threading is properly setup inside DFTK by adding to the script running the DFTK calculation:

    using DFTK
    setup_threading()

    This disables FFT threading and sets the number of BLAS threads to the number of Julia threads.

  2. Run Julia passing the desired number of threads using the flag -t:

    julia -t 8 myscript.jl

For some cases (e.g. a single $k$-point, fewish bands and a large FFT grid) it can be advantageous to add threading inside the FFTs as well. One example is the Caffeine calculation in the above scaling plot. In order to do so just call setup_threading(n_fft=2), which will select two FFT threads. More than two FFT threads is rarely useful.

Advanced threading tweaks

The default threading setup done by setup_threading is to select one FFT thread, and use Threads.nthreads() to determine the number of DFTK and BLAS threads. This section provides some info in case you want to change these defaults.

BLAS threads

All BLAS calls in Julia go through a parallelized OpenBlas or MKL (with MKL.jl. Generally threading in BLAS calls is far from optimal and the default settings can be pretty bad. For example for CPUs with hyper threading enabled, the default number of threads seems to equal the number of virtual cores. Still, BLAS calls typically take second place in terms of the share of runtime they make up (between 10% and 20%). Of note many of these do not take place on matrices of the size of the full FFT grid, but rather only in a subspace (e.g. orthogonalization, Rayleigh-Ritz, ...) such that parallelization is either anyway disabled by the BLAS library or not very effective. To set the number of BLAS threads use

using LinearAlgebra
BLAS.set_num_threads(N)

where N is the number of threads you desire. To check the number of BLAS threads currently used, you can use BLAS.get_num_threads().

DFTK threads

On top of BLAS threading DFTK uses Julia threads in a couple of places to parallelize over $k$-points (density computation) or bands (Hamiltonian application). The number of threads used for these aspects is controlled by default by Threads.nthreads(), which can be set using the flag -t passed to Julia or the environment variable JULIA_NUM_THREADS. Optionally, you can use setup_threading(; n_DFTK) to set an upper bound to the number of threads used by DFTK, regardless of Threads.nthreads() (for instance, if you want to thread several DFTK runs and don't want DFTK's threading to interfere with that).

FFT threads

Since FFT threading is only used in DFTK inside the regions already parallelized by Julia threads, setting FFT threads to something larger than 1 is rarely useful if a sensible number of Julia threads has been chosen. Still, to explicitly set the FFT threads use

using FFTW
FFTW.set_num_threads(N)

where N is the number of threads you desire. By default no FFT threads are used, which is almost always the best choice.