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: 285ms / 54.4% 85.6MiB / 77.5%
Section ncalls time %tot avg alloc %tot avg
────────────────────────────────────────────────────────────────────────────────
self_consistent_field 1 155ms 99.9% 155ms 66.4MiB 100.0% 66.4MiB
LOBPCG 27 66.4ms 42.8% 2.46ms 17.5MiB 26.4% 664KiB
DftHamiltonian... 74 51.1ms 33.0% 691μs 5.86MiB 8.8% 81.1KiB
local+kinetic 489 48.7ms 31.5% 100μs 275KiB 0.4% 576B
nonlocal 74 992μs 0.6% 13.4μs 1.36MiB 2.1% 18.9KiB
rayleigh_ritz 47 4.47ms 2.9% 95.2μs 1.75MiB 2.6% 38.1KiB
ortho! 47 480μs 0.3% 10.2μs 260KiB 0.4% 5.53KiB
ortho! X vs Y 67 4.24ms 2.7% 63.2μs 1.93MiB 2.9% 29.6KiB
ortho! 133 2.07ms 1.3% 15.6μs 1.04MiB 1.6% 7.99KiB
preconditioning 74 989μs 0.6% 13.4μs 274KiB 0.4% 3.71KiB
Update residuals 74 669μs 0.4% 9.04μs 1.13MiB 1.7% 15.6KiB
ortho! 27 474μs 0.3% 17.6μs 166KiB 0.2% 6.16KiB
compute_density 9 57.3ms 37.0% 6.37ms 6.26MiB 9.4% 712KiB
symmetrize_ρ 9 50.7ms 32.7% 5.63ms 5.01MiB 7.5% 570KiB
energy_hamiltonian 19 26.2ms 16.9% 1.38ms 30.5MiB 45.9% 1.61MiB
ene_ops 19 24.0ms 15.5% 1.26ms 19.9MiB 30.0% 1.05MiB
ene_ops: xc 19 20.2ms 13.0% 1.06ms 8.92MiB 13.4% 481KiB
ene_ops: har... 19 2.38ms 1.5% 125μs 8.56MiB 12.9% 462KiB
ene_ops: non... 19 482μs 0.3% 25.4μs 291KiB 0.4% 15.3KiB
ene_ops: kin... 19 364μs 0.2% 19.1μs 180KiB 0.3% 9.46KiB
ene_ops: local 19 250μs 0.2% 13.2μs 1.74MiB 2.6% 93.8KiB
ortho_qr 3 125μs 0.1% 41.7μs 102KiB 0.1% 33.9KiB
χ0Mixing 9 34.4μs 0.0% 3.82μs 37.0KiB 0.1% 4.11KiB
enforce_real! 1 80.0μs 0.1% 80.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 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.
First disable all threading inside DFTK, by adding the following to your script running the DFTK calculation:
using DFTK disable_threading()
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 usempiexecjl
to automatically select thempiexec
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.
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:
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.
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 the same number of BLAS and Julia 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
Int(ccall((BLAS.@blasfunc(openblas_get_num_threads), BLAS.libblas), Cint, ()))
or (from Julia 1.6) simply BLAS.get_num_threads()
.
Julia threads
On top of BLAS threading DFTK uses Julia threads (Thread.@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 the flag -t
passed to Julia or the environment variable JULIA_NUM_THREADS
. To check the number of Julia threads use Threads.nthreads()
.
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.