Important documents and further reading
CI and automatic testing
With every commit and also as nightly checks, continuous integration (CI) tests are run. There are several CI pipelines covering different aspects which are defined in Build -> Pipelines of the project's navigation bar. The typical test is a regression test that runs the code for a minimal example (should run very quickly, not accurately) of a certain feature and checks if the output agrees to some reference.
The tests are defined in the tests
folder where the different stages (build, run and verify steps for different codes) are defined in yaml files collected in tests/gitlab-ci
. A typical step looks like this (annotated for readability):
build_kkrhost:gnu:debug: # name of the job
stage: build_kkrhost # stage where it runs (see pipeline definition)
needs: ["test:intel"] # dependency for this job (here requires to have checked if intel compiler works)
tags:
- docker-executor # required to run the docker image where this is based upon
image: iffregistry.fz-juelich.de/docker-images/centos7-intel-compilers/extended_gcc8:latest # definition of docker image that has intel compilers installed
allow_failure: false
script: # definition of steps that are done for compilation
# first make sure correct cmake etc is used
- source /etc/profile
- mkdir build_gfortran_debug && cd build_gfortran_debug
- /opt/gcc8/bin/gfortran --version
- FC=/opt/gcc8/bin/gfortran cmake -DENABLE_MPI=OFF -DENABLE_OMP=OFF -DCMAKE_BUILD_TYPE=Debug -DCOMPILE_KKRHOST=ON -DCOMPILE_KKRIMP=OFF -DCOMPILE_KKRSUSC=OFF -DCOMPILE_PKKPRIME=OFF -DCOMPILE_VORONOI=OFF -DCOMPILE_RHOQ=OFF ..
- make -j4 VERBOSE=1
only: # definition for which pipelines this test is run
- BdG # branch name
- schedules # scheduled tests (nightly builds)
- triggers # can be triggered by webhooks
- web # can be started from the website
- master # branch name
- develop # branch name
A typical regression test contains the minimal input (e.g., for KKRhost qdos test: inputcard
, shapefun
, potential
, qvec.dat
) and reference data (qdos.AAA.S.dat
output files). After compiling the executable, the CI system then runs the tests and verifies the outcome using simple Python scripts and pytest.
For the qdos test, this is defined in tests/KKRhost/tools/test_verify_qdos.py
.
An overview of the code's test coverage is automatically computed by running the tests with the appropriate compiler option.
Zenodo releases
Code releases are published on Zenodo providing a citable entry: doi: 10.5281/zenodo.7284738.
A new release is automatically triggered when a release is created on the github mirror of JuKKR (the github clone mainly exists for better visibility of the code).
...
in the source code?
Where do I find Here is a loose collection of JuKKR features and where they are actually implemented. This list is not comprehensive but indicates a possible starting point for future development.
Overview of source code
JuKKR collects a number of programs that share common source code (e.g., the Chebychev radial solver for the single site part). The following list gives an overview of the source code structure
-
source/common/
: common code shared over multiple programs -
source/voronoi/
: Voronoi code for structure setup and shape function generation (base code for KKRhost) -
source/KKRhost/
: host code for 2D and 3D periodic crystals (base code for KKRimp, KKRsusc, PKKprime) -
source/KKRimp/
: impurity embedding code into a periodic crystal (base code for KKRsusc) -
source/KKRsusc/
: susceptibility code (TD-DFT) -
source/KKRnano/
: linear scaling code for extremely large unit cells -
source/PKKprime/
: Fermi surface and Boltzmann transport code
Special features
Bogoliubov-de Gennes (BdG) for superconductivity
The Bogoliubov-de Gennes implementation is found mainly on the BdG
branch of JuKKR. Since the BdG implementation requires many changes in matrix sizes, these changes appear throughout the code but they usually are triggered with the kBdG
parameter which changes from kBdG=0
to kBdG=1
when the BdG formalism is used. This then translates into matrix dimension that change to lmmax*(1+kBdG)
taking this into account.
The main changes to the radial solver are published in the KS-BdG method paper.
Useful might also be the following notes on the KKRhost code detailing the KKR multiple scattering formalism with implementation details and important changes introduced by the BdG formalism.
Broyden mixing
Chebychev radial solver
common/radial_solver_Chebychev/rllsll.F90
- see also chapter 5 of PhD thesis by David Bauer
- and also Rudolf Zeller, Frontiers in Physics, 12, (2024). doi: 10.3389/fphy.2024.1393130
Coherent potential approximation (CPA)
Useful information on the CPA algorithm and its implementation in JuKKR including reference points (names of source files found in source/KKRhost
) to the JuKKR specific implementation are summarized in these notes By Manuel dos Santos Dias.
Energy integration contour
J_{ij}, \vec{D}_{ij}
exchange coupling parameters
Extraction of - see the corresponding wiki page
-
source/KKRhost/jijhelp.f90
: calculation of\Delta t
matrices -
source/KKRhost/tbxccpljijdij.F90
: calculation ofJ_{ij}, \vec{D}_{ij}
based on Liechtenstein formula
Inputs for PKKprime code
Files necessary for Fermi surface generation:
-
source/KKRhost/write_tbkkr_files.f90
Additional files needed for Boltzmann transport calculations: source/KKRhost/operators_for_FScode.F90
source/KKRhost/normcoeff_SO.F90
MPI parallelization
Noncollinear magnetism: rotation between local and global spin frame
source/common/rotatespinframe.f90
- see also chapter 6.3 of PhD thesis by David Bauer
inputcard
etc.
Parsing of source/KKRhost/rinput13.F90
-
source/common/runoptions.F90
: available run options -
source/common/global_variables.F90
: definition of global constants (array sizes etc.) -
source/common/constants.f90
: definition of physical constants