GPAW + WannierBerri: Altermagnet MnTe

It is assumed that the reader is familiar with the previous tutorial on bcc Fe with GPAW. In this tutorial we will show a special case, the crystal structure is anti-ferromagnetic. MnTe has two Mn atoms per unit cell with opposite directions of magnetic moments. As a result in spin-up states, the Mn1 d-states are occupied, while the Mn2 d-states are empty. In the spin-down channel the situation is reversed. Thus, to wannierise only the valence bands, one may want to introduce different projections for the two spin channels.

[2]:
import numpy as np
from ase import Atoms
from gpaw import GPAW, PW, FermiDirac, MixerSum
from wannierberri.parallel import Parallel, Serial

parallel_env= Parallel(num_cpus=16)
seed = "MnTe"

a = 4.134
c = 6.652
lattice = a * np.array([[1, 0, 0], [-1 / 2, np.sqrt(3) / 2, 0], [0, 0, c / a]])
positions = np.array(
    [
        [0, 0, 0],
        [0, 0, 1 / 2],
        [1 / 3, 2 / 3, 1 / 4],
        [2 / 3, 1 / 3, 3 / 4],
    ]
)
initializing ray with  {'num_cpus': 16}
2025-10-22 21:08:55,827 INFO worker.py:1918 -- Started a local Ray instance. View the dashboard at http://127.0.0.1:8265 

Step 1: Gpaw calculation

loading the GPAW calculators

If you want to skip the GPAW calculation, you can just load the pre-calculated calculators

[15]:
calc_scf = GPAW(f"gpaw/{seed}-scf.gpw")
calc_nscf = GPAW(f"gpaw/{seed}-nscf-irred-664.gpw")
calc_bands = GPAW(f"gpaw/{seed}-bands-path.gpw")

self-consistent calculation

[ ]:
a = 4.134
c = 6.652
lattice = a * np.array([[1, 0, 0], [-1 / 2, np.sqrt(3) / 2, 0], [0, 0, c / a]])
positions = np.array(
    [
        [0, 0, 0],
        [0, 0, 1 / 2],
        [1 / 3, 2 / 3, 1 / 4],
        [2 / 3, 1 / 3, 3 / 4],
    ]
)
atoms = Atoms("Mn2Te2",
                cell=lattice,
                scaled_positions=positions,
                pbc=[1, 1, 1])

m = 4.7
magmoms = np.zeros(4)
magmoms[0] = +m
magmoms[1] = -m
atoms.set_initial_magnetic_moments(magmoms)

calc = GPAW(
    mode=PW(600),
    xc="PBE",
    kpts={"size": [6, 6, 4], "gamma": True},
    convergence={"density": 1e-6},
    mixer=MixerSum(0.25, 8, 100),
    setups={"Mn": ":d,4.0"},
    txt="MnTe_scf_norelax.txt"
)
atoms.calc = calc
atoms.get_potential_energy()
calc.write("{seed}-scf.gpw", mode="all")

Non-self-consistent calculation in the irreducible Brillouin zone

Now we get the list of irreducible k-points (using irrep.SpaceGroup object) and perform a non-self-consistent calculation only on these k-points. Note, that here we use all symmetries of the crystal structure, even those which are broken by the magnetization direction, because they are still symmetries for each spin channel separately (including the time-reversal symmetry).

[ ]:
from irrep.spacegroup import SpaceGroup
sg = SpaceGroup.from_gpaw(calc)
kpoints_irred = sg.get_irreducible_kpoints_grid([6, 6, 4])
calc_nscf = calc.fixed_density(
    kpts=kpoints_irred,
    symmetry={'symmorphic': False},
    nbands=60,
    convergence={'bands': 40},
    txt=f'{seed}-nscf.txt')
calc_nscf.write(f'{seed}-nscf-irred-664.gpw', mode='all')
typat used for spacegroup detection (accounting magmoms): [26]

Compute the dft bandstructure along a high-symmetry path

This is done to compare with the wannierized bandstructure later.

[ ]:
from wannierberri.grid import Path
# recip_lattice = 2*np.pi * np.linalg.inv(lattice).T
kz = 0.35 / (2*np.pi/c)
path = Path(lattice,
            nodes=[
                [2 / 3, -1 / 3, 0],
                [0, 0, 0],
                [-2 / 3, 1 / 3, 0],
                None,
                [-0.5, 0, kz],
                [0, 0, kz],
                [0.5, 0, kz],
            ],
            labels=[r"${\rm K}\leftarrow$",
                    r"$\Gamma$",
                    r"$\rightarrow{\rm K}$",
                    r"$\overline{\rm M}\leftarrow$",
                    r"$\overline{\Gamma}$",
                    r"$\rightarrow\overline{\rm M}$"],
            length=100)   # length [ Ang] ~= 2*pi/dk
kpoints = path.K_list
print (kpoints.shape)
print (kpoints)

calc_bands = calc_scf.fixed_density(
    nbands=30,
    symmetry='off',
    random=True,
    kpts=kpoints,
    txt=f'{seed}-bands-path.txt',
    convergence={'bands': 24})
calc_bands.write(f'{seed}-bands-path.gpw', mode='all')

(61, 3)
[[ 0.66666667 -0.33333333  0.        ]
 [ 0.625      -0.3125      0.        ]
 [ 0.58333333 -0.29166667  0.        ]
 [ 0.54166667 -0.27083333  0.        ]
 [ 0.5        -0.25        0.        ]
 [ 0.45833333 -0.22916667  0.        ]
 [ 0.41666667 -0.20833333  0.        ]
 [ 0.375      -0.1875      0.        ]
 [ 0.33333333 -0.16666667  0.        ]
 [ 0.29166667 -0.14583333  0.        ]
 [ 0.25       -0.125       0.        ]
 [ 0.20833333 -0.10416667  0.        ]
 [ 0.16666667 -0.08333333  0.        ]
 [ 0.125      -0.0625      0.        ]
 [ 0.08333333 -0.04166667  0.        ]
 [ 0.04166667 -0.02083333  0.        ]
 [ 0.          0.          0.        ]
 [-0.04166667  0.02083333  0.        ]
 [-0.08333333  0.04166667  0.        ]
 [-0.125       0.0625      0.        ]
 [-0.16666667  0.08333333  0.        ]
 [-0.20833333  0.10416667  0.        ]
 [-0.25        0.125       0.        ]
 [-0.29166667  0.14583333  0.        ]
 [-0.33333333  0.16666667  0.        ]
 [-0.375       0.1875      0.        ]
 [-0.41666667  0.20833333  0.        ]
 [-0.45833333  0.22916667  0.        ]
 [-0.5         0.25        0.        ]
 [-0.54166667  0.27083333  0.        ]
 [-0.58333333  0.29166667  0.        ]
 [-0.625       0.3125      0.        ]
 [-0.5         0.          0.37054454]
 [-0.46428571  0.          0.37054454]
 [-0.42857143  0.          0.37054454]
 [-0.39285714  0.          0.37054454]
 [-0.35714286  0.          0.37054454]
 [-0.32142857  0.          0.37054454]
 [-0.28571429  0.          0.37054454]
 [-0.25        0.          0.37054454]
 [-0.21428571  0.          0.37054454]
 [-0.17857143  0.          0.37054454]
 [-0.14285714  0.          0.37054454]
 [-0.10714286  0.          0.37054454]
 [-0.07142857  0.          0.37054454]
 [-0.03571429  0.          0.37054454]
 [ 0.          0.          0.37054454]
 [ 0.03571429  0.          0.37054454]
 [ 0.07142857  0.          0.37054454]
 [ 0.10714286  0.          0.37054454]
 [ 0.14285714  0.          0.37054454]
 [ 0.17857143  0.          0.37054454]
 [ 0.21428571  0.          0.37054454]
 [ 0.25        0.          0.37054454]
 [ 0.28571429  0.          0.37054454]
 [ 0.32142857  0.          0.37054454]
 [ 0.35714286  0.          0.37054454]
 [ 0.39285714  0.          0.37054454]
 [ 0.42857143  0.          0.37054454]
 [ 0.46428571  0.          0.37054454]
 [ 0.5         0.          0.37054454]]
[14]:
bs_dft = calc_bands.band_structure()
bs_dft.plot(show=True, emin=-8,emax=15.0)
../../../_images/tutorials_8_GPAW_3.MnTe_mnte_10_0.png
[14]:
<Axes: ylabel='energies [eV]'>

Step 2: Wannierization

In this case we use sp3d2 and t2g projections on each iron atom. One can also use the conventional s,p,d orbitals.

[16]:
from irrep.spacegroup import SpaceGroup
from wannierberri.symmetry.projections import Projection, ProjectionsSet

sg = SpaceGroup.from_gpaw(calc_nscf)

positions_Mn = [[0, 0, 0],
                [0, 0, 1 / 2]]

positions_Te = [[1 / 3, 2 / 3, 1 / 4],
                [2 / 3, 1 / 3, 3 / 4]]

proj_Mn1_d = Projection(position_num=positions_Mn[0], orbital='d', spacegroup=sg)
proj_Mn2_d = Projection(position_num=positions_Mn[1], orbital='d', spacegroup=sg)

proj_Te_sp2 = Projection(position_num=positions_Te, orbital='sp2', spacegroup=sg, xaxis=[0, -1, 0], rotate_basis=True)
proj_Te_pz = Projection(position_num=positions_Te, orbital='pz', spacegroup=sg)

proj_set_up = ProjectionsSet([proj_Mn1_d, proj_Te_sp2, proj_Te_pz])
proj_set_down = ProjectionsSet([proj_Mn2_d, proj_Te_sp2, proj_Te_pz])


typat used for spacegroup detection (accounting magmoms): [25003, 25001, 52002, 52002]

create the “w90 files”

Here we are NOT using the gpaw-wannier90 interface, and actually not creating the w90 files, but directly access the wavefunctions from the GPAW calculation, use symmetry operations from irrep, and create the necessary data to be used with wannierberri. Those files still retain the same naming convention as the w90 files for consistency, but htey are binary npz files, which are convenient to work with numpy

After creating those data, it is written to the disk, so that one can later load it. But be careful - if you change anything (e.g. the projections) - do not forget to remove those files, otherwise you will load the old data! For that reason, here we set ‘read_npz_list=[]’ to force the creation of new data.

[17]:
from wannierberri.w90files.w90data_soc import Wannier90dataSOC
w90data = Wannier90dataSOC.from_gpaw(
    seedname="wannier_soc",
    calculator=calc_nscf,
    projections_up=proj_set_up,
    projections_down=proj_set_down,
    mp_grid=(6, 6, 4),
    # read_npz_list=[],
    spacegroup=sg,
)

w90data.select_bands(win_min=-10,
                     win_max=50)
finding num points from 3 projections
got irreducible=None, mp_grid=(6, 6, 4), seedname=wannier_soc-spin-0, files=['mmn', 'eig', 'amn', 'symmetrizer'], read_npz_list=None, write_npz_list=None, projections=ProjectionsSet with 13 Wannier functions and 0 free variables
Projection 0, 0, 0:['d'] with 5 Wannier functions on 1 points (5 per site)
Projection 0.3333333333333333, 0.6666666666666666, 0.25:['sp2'] with 6 Wannier functions on 2 points (3 per site)
Projection 0.3333333333333333, 0.6666666666666666, 0.25:['pz'] with 2 Wannier functions on 2 points (1 per site), unk_grid=None, normalize=True
kpt_latt_grid=[[0.         0.         0.        ]
 [0.         0.         0.25      ]
 [0.         0.         0.5       ]
 [0.         0.16666667 0.        ]
 [0.         0.16666667 0.25      ]
 [0.         0.16666667 0.5       ]
 [0.         0.16666667 0.75      ]
 [0.         0.33333333 0.        ]
 [0.         0.33333333 0.25      ]
 [0.         0.33333333 0.5       ]
 [0.         0.33333333 0.75      ]
 [0.         0.5        0.        ]
 [0.         0.5        0.25      ]
 [0.         0.5        0.5       ]
 [0.16666667 0.16666667 0.        ]
 [0.16666667 0.16666667 0.25      ]
 [0.16666667 0.16666667 0.5       ]
 [0.16666667 0.33333333 0.        ]
 [0.16666667 0.33333333 0.25      ]
 [0.16666667 0.33333333 0.5       ]
 [0.16666667 0.33333333 0.75      ]
 [0.33333333 0.33333333 0.        ]
 [0.33333333 0.33333333 0.25      ]
 [0.33333333 0.33333333 0.5       ]]
mpgrid = [6 6 4], 24
detected grid=(np.int64(6), np.int64(6), np.int64(4)), selected_kpoints=[ 0  1  2  3  4  5  6  7  8  9 10 11 12 13 14 15 16 17 18 19 20 21 22 23]
self.irreducible=True
mpgrid = [6 6 4], 24
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/w90files/w90data.py:151: UserWarning: detected grid (np.int64(6), np.int64(6), np.int64(4)) od 144 kpoints, but only 24 kpoints are available.assuming that only irreducible kpoints are needed.
  warnings.warn(f"detected grid {grid} od {np.prod(grid)} kpoints, "
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/w90files/w90data.py:185: UserWarning: Failed to read symmetrizer from wannier_soc-spin-0.symmetrizer.npz: [Errno 2] No such file or directory: 'wannier_soc-spin-0.symmetrizer.npz'
  warnings.warn(f"Failed to read symmetrizer from {fname}: {e}")
orbitals = ['d']
orbitals = ['sp2']
orbitals = ['pz']
calculating Wannier functions for d at [[0 0 0]]
calculating Wannier functions for sp2 at [[0.33333333 0.66666667 0.25      ]
 [0.66666667 0.33333333 0.75      ]]
calculating Wannier functions for pz at [[0.33333333 0.66666667 0.25      ]
 [0.66666667 0.33333333 0.75      ]]
D.shape [(24, 24, 5, 5), (24, 24, 6, 6), (24, 24, 2, 2)]
num_wann 13
D_wann_block_indices [[ 0  5]
 [ 5 11]
 [11 13]]
saving to wannier_soc-spin-0.symmetrizer.npz :
saving to wannier_soc-spin-0.eig.npz :
Creating amn. Using projections_set
ProjectionsSet with 13 Wannier functions and 0 free variables
Projection 0, 0, 0:['d'] with 5 Wannier functions on 1 points (5 per site)
Projection 0.3333333333333333, 0.6666666666666666, 0.25:['sp2'] with 6 Wannier functions on 2 points (3 per site)
Projection 0.3333333333333333, 0.6666666666666666, 0.25:['pz'] with 2 Wannier functions on 2 points (1 per site)
saving to wannier_soc-spin-0.amn.npz :
mpgrid = [6 6 4], 144
NK= 144, selected_kpoints = [ 0  1  2  3  4  5  6  7  8  9 10 11 12 13 14 15 16 17 18 19 20 21 22 23], kptirr = [ 0  1  2  3  4  5  6  7  8  9 10 11 12 13 14 15 16 17 18 19 20 21 22 23]
Shells found with weights [8.96674318 3.89604278] and tolerance 3.147886549272457e-16
saving to wannier_soc-spin-0.mmn.npz :
finding num points from 3 projections
got irreducible=None, mp_grid=(6, 6, 4), seedname=wannier_soc-spin-1, files=['mmn', 'eig', 'amn', 'symmetrizer'], read_npz_list=None, write_npz_list=None, projections=ProjectionsSet with 13 Wannier functions and 0 free variables
Projection 0.0, 0.0, 0.5:['d'] with 5 Wannier functions on 1 points (5 per site)
Projection 0.3333333333333333, 0.6666666666666666, 0.25:['sp2'] with 6 Wannier functions on 2 points (3 per site)
Projection 0.3333333333333333, 0.6666666666666666, 0.25:['pz'] with 2 Wannier functions on 2 points (1 per site), unk_grid=None, normalize=True
kpt_latt_grid=[[0.         0.         0.        ]
 [0.         0.         0.25      ]
 [0.         0.         0.5       ]
 [0.         0.16666667 0.        ]
 [0.         0.16666667 0.25      ]
 [0.         0.16666667 0.5       ]
 [0.         0.16666667 0.75      ]
 [0.         0.33333333 0.        ]
 [0.         0.33333333 0.25      ]
 [0.         0.33333333 0.5       ]
 [0.         0.33333333 0.75      ]
 [0.         0.5        0.        ]
 [0.         0.5        0.25      ]
 [0.         0.5        0.5       ]
 [0.16666667 0.16666667 0.        ]
 [0.16666667 0.16666667 0.25      ]
 [0.16666667 0.16666667 0.5       ]
 [0.16666667 0.33333333 0.        ]
 [0.16666667 0.33333333 0.25      ]
 [0.16666667 0.33333333 0.5       ]
 [0.16666667 0.33333333 0.75      ]
 [0.33333333 0.33333333 0.        ]
 [0.33333333 0.33333333 0.25      ]
 [0.33333333 0.33333333 0.5       ]]
mpgrid = [6 6 4], 24
detected grid=(np.int64(6), np.int64(6), np.int64(4)), selected_kpoints=[ 0  1  2  3  4  5  6  7  8  9 10 11 12 13 14 15 16 17 18 19 20 21 22 23]
self.irreducible=True
mpgrid = [6 6 4], 24
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/w90files/w90data.py:185: UserWarning: Failed to read symmetrizer from wannier_soc-spin-1.symmetrizer.npz: [Errno 2] No such file or directory: 'wannier_soc-spin-1.symmetrizer.npz'
  warnings.warn(f"Failed to read symmetrizer from {fname}: {e}")
orbitals = ['d']
orbitals = ['sp2']
orbitals = ['pz']
calculating Wannier functions for d at [[0.  0.  0.5]]
calculating Wannier functions for sp2 at [[0.33333333 0.66666667 0.25      ]
 [0.66666667 0.33333333 0.75      ]]
calculating Wannier functions for pz at [[0.33333333 0.66666667 0.25      ]
 [0.66666667 0.33333333 0.75      ]]
D.shape [(24, 24, 5, 5), (24, 24, 6, 6), (24, 24, 2, 2)]
num_wann 13
D_wann_block_indices [[ 0  5]
 [ 5 11]
 [11 13]]
saving to wannier_soc-spin-1.symmetrizer.npz :
saving to wannier_soc-spin-1.eig.npz :
Creating amn. Using projections_set
ProjectionsSet with 13 Wannier functions and 0 free variables
Projection 0.0, 0.0, 0.5:['d'] with 5 Wannier functions on 1 points (5 per site)
Projection 0.3333333333333333, 0.6666666666666666, 0.25:['sp2'] with 6 Wannier functions on 2 points (3 per site)
Projection 0.3333333333333333, 0.6666666666666666, 0.25:['pz'] with 2 Wannier functions on 2 points (1 per site)
saving to wannier_soc-spin-1.amn.npz :
mpgrid = [6 6 4], 144
NK= 144, selected_kpoints = [ 0  1  2  3  4  5  6  7  8  9 10 11 12 13 14 15 16 17 18 19 20 21 22 23], kptirr = [ 0  1  2  3  4  5  6  7  8  9 10 11 12 13 14 15 16 17 18 19 20 21 22 23]
Shells found with weights [8.96674318 3.89604278] and tolerance 3.147886549272457e-16
saving to wannier_soc-spin-1.mmn.npz :
number of bands = 60
saving to wannier_soc.soc.npz :
selected_bands = [ 8  9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31
 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55
 56 57 58 59]
selected_bands = [ 8  9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31
 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55
 56 57 58 59]

Do the wannierization

at this step we also specify the frozen window, number of iterations, and whether to use site symmetry and localization. (when starting from irreducible BZ, site-symmetry is required anyway)

We also may use the outer window via w90data.select_bands if needed

In this case, both spin channels will be wannierised separately, starting from the same projections.There is also an option to specify different projections for different spin channels if needed. (See MnTe tutorial)

[19]:
w90data.wannierise(
    froz_min=-10,
    froz_max=7,
    num_iter=300,
    print_progress_every=50,
    sitesym=True,
    localise=True,
)


Symmetrizer_Uirr initialized for ikirr=0, kpt=0, [0. 0. 0.] with 24 symmetries, max error in included blocks: 4.775301491256433e-07 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 2.8627669436984223e-06, 9.576384343585948e-07
Symmetrizer_Uirr initialized for ikirr=1, kpt=1, [0.   0.   0.25] with 12 symmetries, max error in included blocks: 3.381164304429173e-06 ; excluded bands are [48 49] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=2, kpt=2, [0.  0.  0.5] with 24 symmetries, max error in included blocks: 1.974197330420772e-06 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=3, kpt=3, [0.         0.16666667 0.        ] with 4 symmetries, max error in included blocks: 2.799612726714073e-07 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 33 [47:48] is 3.381164304429173e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 31 [46:48] is 1.974197330420772e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
Warning: symmetrization did not converge in 100 iterations, final changes 0.00012340121863906617, 4.1580270690599557e-07
Symmetrizer_Uirr initialized for ikirr=4, kpt=4, [0.         0.16666667 0.25      ] with 4 symmetries, max error in included blocks: 3.7372458106401275e-07 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 0.00015796458368711307, 4.450804368802042e-05
Symmetrizer_Uirr initialized for ikirr=5, kpt=5, [0.         0.16666667 0.5       ] with 4 symmetries, max error in included blocks: 9.306464534587413e-07 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 6.660103712958651e-05, 6.513289194598715e-07
Symmetrizer_Uirr initialized for ikirr=6, kpt=6, [0.         0.16666667 0.75      ] with 4 symmetries, max error in included blocks: 8.359040398240567e-07 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 0.00010263731306014524, 2.6830872368936694e-05
Symmetrizer_Uirr initialized for ikirr=7, kpt=7, [0.         0.33333333 0.        ] with 4 symmetries, max error in included blocks: 2.9180980241263033e-08 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 1.4601151838923914e-05, 2.892621532286757e-06
Symmetrizer_Uirr initialized for ikirr=8, kpt=8, [0.         0.33333333 0.25      ] with 4 symmetries, max error in included blocks: 5.038406993695107e-07 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 5.773736105965039e-05, 1.0816771404676023e-05
Symmetrizer_Uirr initialized for ikirr=9, kpt=9, [0.         0.33333333 0.5       ] with 4 symmetries, max error in included blocks: 4.656445783754465e-07 ; excluded bands are [48 49 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 2.752614722251489e-06, 7.532731462939769e-07
Symmetrizer_Uirr initialized for ikirr=10, kpt=10, [0.         0.33333333 0.75      ] with 4 symmetries, max error in included blocks: 7.920167944688993e-07 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 0.0010894560510402904, 0.00023541711487988964
Symmetrizer_Uirr initialized for ikirr=11, kpt=11, [0.  0.5 0. ] with 8 symmetries, max error in included blocks: 3.57990449657817e-06 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=12, kpt=12, [0.   0.5  0.25] with 4 symmetries, max error in included blocks: 6.099941322125381e-07 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 47 [47:48] is 3.57990449657817e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
Warning: symmetrization did not converge in 100 iterations, final changes 2.6248060891252107e-05, 3.1836627036827746e-06
Symmetrizer_Uirr initialized for ikirr=13, kpt=13, [0.  0.5 0.5] with 8 symmetries, max error in included blocks: 5.976570863867099e-08 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 1.1292944388973457e-06, 3.3077916259809455e-07
Symmetrizer_Uirr initialized for ikirr=14, kpt=14, [0.16666667 0.16666667 0.        ] with 4 symmetries, max error in included blocks: 4.180302175212303e-07 ; excluded bands are [50 51] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=15, kpt=15, [0.16666667 0.16666667 0.25      ] with 2 symmetries, max error in included blocks: 1.1739862079479962e-16 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 3.96769961262693e-06, 1.0576300232120623e-06
Symmetrizer_Uirr initialized for ikirr=16, kpt=16, [0.16666667 0.16666667 0.5       ] with 4 symmetries, max error in included blocks: 6.688390301510412e-07 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=17, kpt=17, [0.16666667 0.33333333 0.        ] with 2 symmetries, max error in included blocks: 1.1857187100668868e-16 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=18, kpt=18, [0.16666667 0.33333333 0.25      ] with 2 symmetries, max error in included blocks: 1.0007415106216804e-16 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=19, kpt=19, [0.16666667 0.33333333 0.5       ] with 2 symmetries, max error in included blocks: 1.1188630228279524e-16 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=20, kpt=20, [0.16666667 0.33333333 0.75      ] with 2 symmetries, max error in included blocks: 1.2412670766236366e-16 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 0.0012291523968025627, 0.00013047366390576224
Symmetrizer_Uirr initialized for ikirr=21, kpt=21, [0.33333333 0.33333333 0.        ] with 12 symmetries, max error in included blocks: 0.002192262335905529 ; excluded bands are [49] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=22, kpt=22, [0.33333333 0.33333333 0.25      ] with 6 symmetries, max error in included blocks: 1.5244972126722465e-05 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=23, kpt=23, [0.33333333 0.33333333 0.5       ] with 12 symmetries, max error in included blocks: 7.225655578238547e-07 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 9 [14:16] is 1.7675104934565142e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 28 [43:45] is 1.7130895244994746e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 29 [45:46] is 0.0007906085147223288, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 30 [46:48] is 0.002192262335905529, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 11 [17:19] is 1.5244972126722465e-05, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
####################################################################################################
starting WFs
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    0.000000000000   |     0.453612612686
  0.000000000000    0.000000000000    0.000000000000   |     0.488675433049
  0.000000000000    0.000000000000    0.000000000000   |     0.488675433049
  0.000000000000    0.000000000000    0.000000000000   |     0.470020830251
  0.000000000000    0.000000000000    0.000000000000   |     0.470020830251
  0.701560854875    2.791812361245    1.687986300531   |     2.353922087855
 -0.701560854875    2.791812361245    1.687986300531   |     2.353922087855
  0.000000000000    1.576673316000    1.687986300531   |     2.353922087855
  1.365439145125    0.788336658000    4.964013699469   |     2.353922087855
  2.768560854875    0.788336658000    4.964013699469   |     2.353922087855
  2.067000000000    2.003475703245    4.964013699469   |     2.353922087855
  0.000000000000    2.386766012830    1.715208132269   |     3.700047650893
  2.067000000000    1.193383006415    4.936791867731   |     3.700047650893
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   26.608000000000   |    23.894632968205 <- sum
                                          maximal spread =   3.700047650893
####################################################################################################
####################################################################################################
Iteration 0 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    0.000000000000   |     0.452098928388
  0.000000000000    0.000000000000    0.000000000000   |     0.484514010644
  0.000000000000    0.000000000000    0.000000000000   |     0.484514010644
  0.000000000000    0.000000000000    0.000000000000   |     0.466583796377
  0.000000000000    0.000000000000    0.000000000000   |     0.466583796377
  0.703635688420    2.793010266951    1.689366328431   |     2.323737882370
 -0.703635688420    2.793010266951    1.689366328431   |     2.323737882370
 -0.000000000000    1.574277504588    1.689366328431   |     2.323737882370
  1.363364311580    0.787138752294    4.962633671569   |     2.323737882370
  2.770635688420    0.787138752294    4.962633671569   |     2.323737882370
  2.067000000000    2.005871514657    4.962633671569   |     2.323737882370
 -0.000000000000    2.386766012830    1.716180125324   |     3.680145552428
  2.067000000000    1.193383006415    4.935819874676   |     3.680145552428
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   26.608000000000   |    23.657012941506 <- sum
                                          maximal spread =   3.680145552428
standard deviation = 0.0
####################################################################################################
####################################################################################################
Iteration 50 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    0.000000000000   |     0.453703462911
  0.000000000000    0.000000000000    0.000000000000   |     0.476621333689
  0.000000000000    0.000000000000    0.000000000000   |     0.476621333689
  0.000000000000    0.000000000000    0.000000000000   |     0.452960802546
  0.000000000000    0.000000000000    0.000000000000   |     0.452960802546
  0.711601410469    2.797609278720    1.693780673033   |     2.077790598228
 -0.711601410469    2.797609278720    1.693780673033   |     2.077790598228
  0.000000000000    1.565079481050    1.693780673033   |     2.077790598228
  1.355398589531    0.782539740525    4.958219326967   |     2.077790598228
  2.778601410469    0.782539740525    4.958219326967   |     2.077790598228
  2.067000000000    2.015069538195    4.958219326967   |     2.077790598228
 -0.000000000000    2.386766012830    1.732070362801   |     3.551858103671
  2.067000000000    1.193383006415    4.919929637199   |     3.551858103671
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   26.608000000000   |    21.883327532092 <- sum
                                          maximal spread =   3.551858103671
standard deviation = 0.00022616806717684447
####################################################################################################
####################################################################################################
Iteration 100 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    0.000000000000   |     0.453970051338
  0.000000000000    0.000000000000    0.000000000000   |     0.479520443462
  0.000000000000    0.000000000000    0.000000000000   |     0.479520443462
  0.000000000000    0.000000000000    0.000000000000   |     0.453772941484
  0.000000000000    0.000000000000    0.000000000000   |     0.453772941484
  0.711505084832    2.797553665087    1.688390553505   |     2.077717516040
 -0.711505084832    2.797553665087    1.688390553505   |     2.077717516040
  0.000000000000    1.565190708315    1.688390553505   |     2.077717516040
  1.355494915168    0.782595354157    4.963609446495   |     2.077717516040
  2.778505084832    0.782595354157    4.963609446495   |     2.077717516040
  2.067000000000    2.014958310930    4.963609446495   |     2.077717516040
  0.000000000000    2.386766012830    1.749022753278   |     3.547553828312
  2.067000000000    1.193383006415    4.902977246722   |     3.547553828312
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   26.608000000000   |    21.881969574095 <- sum
                                          maximal spread =   3.547553828312
standard deviation = 0.0003307040471834225
####################################################################################################
####################################################################################################
Iteration 150 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    0.000000000000   |     0.453834367875
  0.000000000000    0.000000000000    0.000000000000   |     0.479840393520
  0.000000000000    0.000000000000    0.000000000000   |     0.479840393520
  0.000000000000    0.000000000000    0.000000000000   |     0.453709568246
  0.000000000000    0.000000000000    0.000000000000   |     0.453709568246
  0.711174909493    2.797363038267    1.680421471235   |     2.079612749398
 -0.711174909493    2.797363038267    1.680421471235   |     2.079612749398
 -0.000000000000    1.565571961957    1.680421471235   |     2.079612749398
  1.355825090507    0.782785980978    4.971578528765   |     2.079612749398
  2.778174909493    0.782785980978    4.971578528765   |     2.079612749398
  2.067000000000    2.014577057288    4.971578528765   |     2.079612749398
  0.000000000000    2.386766012830    1.773859521557   |     3.540252224467
  2.067000000000    1.193383006415    4.878140478443   |     3.540252224467
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   26.608000000000   |    21.879115236730 <- sum
                                          maximal spread =   3.540252224467
standard deviation = 0.00048320545104156755
####################################################################################################
####################################################################################################
Iteration 200 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    0.000000000000   |     0.453614303692
  0.000000000000    0.000000000000    0.000000000000   |     0.480143538021
  0.000000000000    0.000000000000    0.000000000000   |     0.480143538021
  0.000000000000    0.000000000000    0.000000000000   |     0.453579809677
  0.000000000000    0.000000000000    0.000000000000   |     0.453579809677
  0.710411429026    2.796922242614    1.668863081480   |     2.083332884787
 -0.710411429026    2.796922242614    1.668863081480   |     2.083332884787
  0.000000000000    1.566453553263    1.668863081480   |     2.083332884787
  1.356588570974    0.783226776631    4.983136918520   |     2.083332884787
  2.777411429026    0.783226776631    4.983136918520   |     2.083332884787
  2.067000000000    2.013695465982    4.983136918520   |     2.083332884787
 -0.000000000000    2.386766012830    1.809839528675   |     3.525930733890
  2.067000000000    1.193383006415    4.842160471325   |     3.525930733890
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   26.608000000000   |    21.872919775587 <- sum
                                          maximal spread =   3.525930733890
standard deviation = 0.0006936985772773872
####################################################################################################
####################################################################################################
Iteration 250 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    0.000000000000   |     0.453308216807
  0.000000000000    0.000000000000    0.000000000000   |     0.480548726224
  0.000000000000    0.000000000000    0.000000000000   |     0.480548726224
  0.000000000000    0.000000000000    0.000000000000   |     0.453415886164
  0.000000000000    0.000000000000    0.000000000000   |     0.453415886164
  0.708751098008    2.795963650053    1.652536025627   |     2.090472317163
 -0.708751098008    2.795963650053    1.652536025627   |     2.090472317163
 -0.000000000000    1.568370738384    1.652536025627   |     2.090472317163
  1.358248901992    0.784185369192    4.999463974373   |     2.090472317163
  2.775751098008    0.784185369192    4.999463974373   |     2.090472317163
  2.067000000000    2.011778280861    4.999463974373   |     2.090472317163
 -0.000000000000    2.386766012830    1.860577476769   |     3.498046730560
  2.067000000000    1.193383006415    4.791422523231   |     3.498046730560
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   26.608000000000   |    21.860164805680 <- sum
                                          maximal spread =   3.498046730560
standard deviation = 0.0009602853117751312
####################################################################################################
####################################################################################################
Final state (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    0.000000000000   |     0.452917212346
  0.000000000000    0.000000000000    0.000000000000   |     0.481046823833
  0.000000000000    0.000000000000    0.000000000000   |     0.481046823833
  0.000000000000    0.000000000000    0.000000000000   |     0.453241413466
  0.000000000000    0.000000000000    0.000000000000   |     0.453241413466
  0.705488709520    2.794080109182    1.631128511626   |     2.103344390200
 -0.705488709520    2.794080109182    1.631128511626   |     2.103344390200
  0.000000000000    1.572137820126    1.631128511626   |     2.103344390200
  1.361511290480    0.786068910063    5.020871488374   |     2.103344390200
  2.772488709520    0.786068910063    5.020871488374   |     2.103344390200
  2.067000000000    2.008011199119    5.020871488374   |     2.103344390200
  0.000000000000    2.386766012830    1.926898514964   |     3.447634425104
  2.067000000000    1.193383006415    4.725101485036   |     3.447634425104
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   26.608000000000   |    21.836828878352 <- sum
                                          maximal spread =   3.447634425104
standard deviation = 0.0012338873736273156
####################################################################################################
time for creating wannierizer 2.207609176635742
time for iterations 11.54856562614441
time for updating 8.583920240402222
total time for wannierization 15.094946384429932
saving to wannier_soc-spin-0.chk :
Symmetrizer_Uirr initialized for ikirr=0, kpt=0, [0. 0. 0.] with 24 symmetries, max error in included blocks: 4.214831298604167e-07 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=1, kpt=1, [0.   0.   0.25] with 12 symmetries, max error in included blocks: 1.3130784851360061e-06 ; excluded bands are [48 50 51] out of 52 total bands (accuracy threshold 1e-06)
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 33 [47:48] is 1.3130784851360061e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
Symmetrizer_Uirr initialized for ikirr=2, kpt=2, [0.  0.  0.5] with 24 symmetries, max error in included blocks: 5.736129865047892e-07 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=3, kpt=3, [0.         0.16666667 0.        ] with 4 symmetries, max error in included blocks: 9.579951910506201e-09 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 3.7639962856174094e-05, 7.28036024830578e-06
Symmetrizer_Uirr initialized for ikirr=4, kpt=4, [0.         0.16666667 0.25      ] with 4 symmetries, max error in included blocks: 4.305908219256872e-08 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 3.652439230064339e-05, 7.1105190068432514e-06
Symmetrizer_Uirr initialized for ikirr=5, kpt=5, [0.         0.16666667 0.5       ] with 4 symmetries, max error in included blocks: 1.842572016762308e-07 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 0.00020026021939457514, 6.374590616106042e-05
Symmetrizer_Uirr initialized for ikirr=6, kpt=6, [0.         0.16666667 0.75      ] with 4 symmetries, max error in included blocks: 4.115514476001744e-07 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 0.00011959051721699526, 1.8924614759873406e-05
Symmetrizer_Uirr initialized for ikirr=7, kpt=7, [0.         0.33333333 0.        ] with 4 symmetries, max error in included blocks: 1.0258829566124216e-07 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 0.0015204181769065501, 0.0003162760082778553
Symmetrizer_Uirr initialized for ikirr=8, kpt=8, [0.         0.33333333 0.25      ] with 4 symmetries, max error in included blocks: 1.7371188884375128e-06 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 0.00088169951917363, 0.0001365773911119321
Symmetrizer_Uirr initialized for ikirr=9, kpt=9, [0.         0.33333333 0.5       ] with 4 symmetries, max error in included blocks: 7.449124970218909e-08 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 45 [46:47] is 1.7371188884375128e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
Warning: symmetrization did not converge in 100 iterations, final changes 0.0007712566764268882, 0.00011800879966863224
Symmetrizer_Uirr initialized for ikirr=10, kpt=10, [0.         0.33333333 0.75      ] with 4 symmetries, max error in included blocks: 3.799296949479876e-06 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 46 [46:47] is 3.799296949479876e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
Warning: symmetrization did not converge in 100 iterations, final changes 0.00018680081651648053, 3.191453704736904e-05
Symmetrizer_Uirr initialized for ikirr=11, kpt=11, [0.  0.5 0. ] with 8 symmetries, max error in included blocks: 1.191229412573126e-06 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=12, kpt=12, [0.   0.5  0.25] with 4 symmetries, max error in included blocks: 1.3824749964370088e-07 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 47 [47:48] is 1.191229412573126e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
Warning: symmetrization did not converge in 100 iterations, final changes 0.0005362438601545182, 0.0002025021029274615
Symmetrizer_Uirr initialized for ikirr=13, kpt=13, [0.  0.5 0.5] with 8 symmetries, max error in included blocks: 7.302686515522201e-06 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=14, kpt=14, [0.16666667 0.16666667 0.        ] with 4 symmetries, max error in included blocks: 6.279557009874495e-07 ; excluded bands are [49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=15, kpt=15, [0.16666667 0.16666667 0.25      ] with 2 symmetries, max error in included blocks: 1.3947004136484323e-16 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=16, kpt=16, [0.16666667 0.16666667 0.5       ] with 4 symmetries, max error in included blocks: 1.1403355381199705e-08 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=17, kpt=17, [0.16666667 0.33333333 0.        ] with 2 symmetries, max error in included blocks: 1.1102230246251565e-16 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 46 [46:47] is 1.018895160721902e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 47 [47:48] is 7.302686515522201e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
Symmetrizer_Uirr initialized for ikirr=18, kpt=18, [0.16666667 0.33333333 0.25      ] with 2 symmetries, max error in included blocks: 1.1857187100668868e-16 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=19, kpt=19, [0.16666667 0.33333333 0.5       ] with 2 symmetries, max error in included blocks: 1.2412670766236366e-16 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=20, kpt=20, [0.16666667 0.33333333 0.75      ] with 2 symmetries, max error in included blocks: 1.1775693440128314e-16 ; excluded bands are [] out of 52 total bands (accuracy threshold 1e-06)
Warning: symmetrization did not converge in 100 iterations, final changes 0.0006392619162661863, 0.00017922490530220204
Symmetrizer_Uirr initialized for ikirr=21, kpt=21, [0.33333333 0.33333333 0.        ] with 12 symmetries, max error in included blocks: 0.0031316702788573665 ; excluded bands are [48 49] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=22, kpt=22, [0.33333333 0.33333333 0.25      ] with 6 symmetries, max error in included blocks: 1.8953455214112623e-05 ; excluded bands are [48 49 50 51] out of 52 total bands (accuracy threshold 1e-06)
Symmetrizer_Uirr initialized for ikirr=23, kpt=23, [0.33333333 0.33333333 0.5       ] with 12 symmetries, max error in included blocks: 7.089386334792644e-07 ; excluded bands are [48] out of 52 total bands (accuracy threshold 1e-06)
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 9 [14:16] is 1.4487757124287662e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 28 [43:45] is 2.035405963138839e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 29 [45:46] is 0.0031316702788573665, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 30 [46:48] is 0.003077322062227622, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 11 [17:19] is 1.8953455214112623e-05, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 29 [45:47] is 1.2930043604619416e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sawf_kirr.py:165: UserWarning: Warning: max error in block 30 [47:48] is 1.4447505013695027e-06, exceeding threshold 1e-06, and this is not among the  upper bands(48:52) bands, this may indicate inaccuracy in the input data
  warnings.warn(f"Warning: max error in block {i} [{start}:{end}] is {max_error_in_blocks[i]}, exceeding threshold {accuracy_threshold}, and this is not among the  upper"
####################################################################################################
starting WFs
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    3.326000000000   |     0.453599406127
  0.000000000000    0.000000000000    3.326000000000   |     0.488672291672
  0.000000000000    0.000000000000    3.326000000000   |     0.488672291672
  0.000000000000    0.000000000000    3.326000000000   |     0.470025991047
  0.000000000000    0.000000000000    3.326000000000   |     0.470025991047
  0.701560487452    2.791812149113    1.637967928255   |     2.353924094187
 -0.701560487452    2.791812149113    1.637967928255   |     2.353924094187
 -0.000000000000    1.576673740264    1.637967928255   |     2.353924094187
  1.365439512548    0.788336870132    5.014032071745   |     2.353924094187
  2.768560487452    0.788336870132    5.014032071745   |     2.353924094187
  2.067000000000    2.003475278981    5.014032071745   |     2.353924094187
  0.000000000000    2.386766012830    1.610814863431   |     3.700071607954
  2.067000000000    1.193383006415    5.041185136569   |     3.700071607954
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   43.238000000000   |    23.894683752593 <- sum
                                          maximal spread =   3.700071607954
####################################################################################################
####################################################################################################
Iteration 0 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    3.326000000000   |     0.452085772301
  0.000000000000    0.000000000000    3.326000000000   |     0.484510619823
  0.000000000000    0.000000000000    3.326000000000   |     0.484510619823
  0.000000000000    0.000000000000    3.326000000000   |     0.466589084554
  0.000000000000    0.000000000000    3.326000000000   |     0.466589084554
  0.703635297648    2.793010041338    1.636587650752   |     2.323740824971
 -0.703635297648    2.793010041338    1.636587650752   |     2.323740824971
 -0.000000000000    1.574277955813    1.636587650752   |     2.323740824971
  1.363364702352    0.787138977907    5.015412349248   |     2.323740824971
  2.770635297648    0.787138977907    5.015412349248   |     2.323740824971
  2.067000000000    2.005871063432    5.015412349248   |     2.323740824971
  0.000000000000    2.386766012830    1.609843634937   |     3.680170325704
  2.067000000000    1.193383006415    5.042156365063   |     3.680170325704
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   43.238000000000   |    23.657070782286 <- sum
                                          maximal spread =   3.680170325704
standard deviation = 0.0
####################################################################################################
####################################################################################################
Iteration 50 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    3.326000000000   |     0.453715162072
  0.000000000000    0.000000000000    3.326000000000   |     0.476626841816
  0.000000000000    0.000000000000    3.326000000000   |     0.476626841816
  0.000000000000    0.000000000000    3.326000000000   |     0.452972432584
  0.000000000000    0.000000000000    3.326000000000   |     0.452972432584
  0.711601098998    2.797609098892    1.632183545131   |     2.077802973841
 -0.711601098998    2.797609098892    1.632183545131   |     2.077802973841
  0.000000000000    1.565079840706    1.632183545131   |     2.077802973841
  1.355398901002    0.782539920353    5.019816454869   |     2.077802973841
  2.778601098998    0.782539920353    5.019816454869   |     2.077802973841
  2.067000000000    2.015069178539    5.019816454869   |     2.077802973841
  0.000000000000    2.386766012830    1.594019678354   |     3.551853298394
  2.067000000000    1.193383006415    5.057980321646   |     3.551853298394
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   43.238000000000   |    21.883438150705 <- sum
                                          maximal spread =   3.551853298394
standard deviation = 0.0002249024450025175
####################################################################################################
####################################################################################################
Iteration 100 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    3.326000000000   |     0.453983073552
  0.000000000000    0.000000000000    3.326000000000   |     0.479525940499
  0.000000000000    0.000000000000    3.326000000000   |     0.479525940499
  0.000000000000    0.000000000000    3.326000000000   |     0.453784972109
  0.000000000000    0.000000000000    3.326000000000   |     0.453784972109
  0.711506105169    2.797554254179    1.637543198686   |     2.077722846346
 -0.711506105169    2.797554254179    1.637543198686   |     2.077722846346
  0.000000000000    1.565189530131    1.637543198686   |     2.077722846346
  1.355493894831    0.782594765065    5.014456801314   |     2.077722846346
  2.778506105169    0.782594765065    5.014456801314   |     2.077722846346
  2.067000000000    2.014959489114    5.014456801314   |     2.077722846346
 -0.000000000000    2.386766012830    1.577162264583   |     3.547576918262
  2.067000000000    1.193383006415    5.074837735417   |     3.547576918262
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   43.238000000000   |    21.882095813368 <- sum
                                          maximal spread =   3.547576918262
standard deviation = 0.0003288578056379191
####################################################################################################
####################################################################################################
Iteration 150 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    3.326000000000   |     0.453848095038
  0.000000000000    0.000000000000    3.326000000000   |     0.479845124298
  0.000000000000    0.000000000000    3.326000000000   |     0.479845124298
  0.000000000000    0.000000000000    3.326000000000   |     0.453721787047
  0.000000000000    0.000000000000    3.326000000000   |     0.453721787047
  0.711179205330    2.797365518469    1.645468213372   |     2.079603061873
 -0.711179205330    2.797365518469    1.645468213372   |     2.079603061873
  0.000000000000    1.565567001551    1.645468213372   |     2.079603061873
  1.355820794670    0.782783500775    5.006531786628   |     2.079603061873
  2.778179205330    0.782783500775    5.006531786628   |     2.079603061873
  2.067000000000    2.014582017694    5.006531786628   |     2.079603061873
  0.000000000000    2.386766012830    1.552462592102   |     3.540335579880
  2.067000000000    1.193383006415    5.099537407898   |     3.540335579880
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   43.238000000000   |    21.879271448724 <- sum
                                          maximal spread =   3.540335579880
standard deviation = 0.0004805711926900684
####################################################################################################
####################################################################################################
Iteration 200 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    3.326000000000   |     0.453628912262
  0.000000000000    0.000000000000    3.326000000000   |     0.480147176703
  0.000000000000    0.000000000000    3.326000000000   |     0.480147176703
  0.000000000000    0.000000000000    3.326000000000   |     0.453592269274
  0.000000000000    0.000000000000    3.326000000000   |     0.453592269274
  0.710423429959    2.796929171356    1.656965085743   |     2.083292075177
 -0.710423429959    2.796929171356    1.656965085743   |     2.083292075177
  0.000000000000    1.566439695779    1.656965085743   |     2.083292075177
  1.356576570041    0.783219847889    4.995034914257   |     2.083292075177
  2.777423429959    0.783219847889    4.995034914257   |     2.083292075177
  2.067000000000    2.013709323466    4.995034914257   |     2.083292075177
  0.000000000000    2.386766012830    1.516673521182   |     3.526139809858
  2.067000000000    1.193383006415    5.135326478818   |     3.526139809858
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   43.238000000000   |    21.873139874993 <- sum
                                          maximal spread =   3.526139809858
standard deviation = 0.000690123989941282
####################################################################################################
####################################################################################################
Iteration 250 (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    3.326000000000   |     0.453323859541
  0.000000000000    0.000000000000    3.326000000000   |     0.480551082297
  0.000000000000    0.000000000000    3.326000000000   |     0.480551082297
  0.000000000000    0.000000000000    3.326000000000   |     0.453428493391
  0.000000000000    0.000000000000    3.326000000000   |     0.453428493391
  0.708779504609    2.795980050612    1.673212592858   |     2.090370444988
 -0.708779504609    2.795980050612    1.673212592858   |     2.090370444988
  0.000000000000    1.568337937266    1.673212592858   |     2.090370444988
  1.358220495391    0.784168968633    4.978787407142   |     2.090370444988
  2.775779504609    0.784168968633    4.978787407142   |     2.090370444988
  2.067000000000    2.011811081979    4.978787407142   |     2.090370444988
  0.000000000000    2.386766012830    1.466181268894   |     3.498501430924
  2.067000000000    1.193383006415    5.185818731106   |     3.498501430924
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   43.238000000000   |    21.860508542697 <- sum
                                          maximal spread =   3.498501430924
standard deviation = 0.0009559487757346814
####################################################################################################
####################################################################################################
Final state (from wannierizer)
----------------------------------------------------------------------------------------------------
wannier centers and spreads
----------------------------------------------------------------------------------------------------
  0.000000000000    0.000000000000    3.326000000000   |     0.452933754168
  0.000000000000    0.000000000000    3.326000000000   |     0.481048194924
  0.000000000000    0.000000000000    3.326000000000   |     0.481048194924
  0.000000000000    0.000000000000    3.326000000000   |     0.453253735568
  0.000000000000    0.000000000000    3.326000000000   |     0.453253735568
  0.705546670356    2.794113572886    1.694534999525   |     2.103138483728
 -0.705546670356    2.794113572886    1.694534999525   |     2.103138483728
  0.000000000000    1.572070892718    1.694534999525   |     2.103138483728
  1.361453329644    0.786035446359    4.957465000475   |     2.103138483728
  2.772546670356    0.786035446359    4.957465000475   |     2.103138483728
  2.067000000000    2.008078126527    4.957465000475   |     2.103138483728
  0.000000000000    2.386766012830    1.400120244503   |     3.448498705901
  2.067000000000    1.193383006415    5.251879755497   |     3.448498705901
----------------------------------------------------------------------------------------------------
  8.268000000000   14.320596076979   43.238000000000   |    21.837365929322 <- sum
                                          maximal spread =   3.448498705901
standard deviation = 0.0012298102324473455
####################################################################################################
time for creating wannierizer 2.31028413772583
time for iterations 11.139642238616943
time for updating 8.321329593658447
total time for wannierization 14.737658023834229
saving to wannier_soc-spin-1.chk :

Create the System object to be used in WannierBerri calculations

It is called “System_w90” class, although it now does not use the w90 wannierization, but the class is the same.

Unlike the non-SOC systems, here we also specify the direction of the magnetization. We can later easily chanbge this direction without redoing the wannierization.

[20]:
from wannierberri.system.system_soc import SystemSOC
theta = 90
phi = 90

system_soc = SystemSOC.from_wannier90data_soc(w90data=w90data, berry=True, silent=False)
system_soc.set_soc_axis(theta=theta, phi=phi, alpha_soc=1.0, units="degrees")
system_soc.save_npz("system_soc")

setting Rvec
setting AA..
setting AA - OK
Real-space lattice:
 [[ 4.134       0.          0.        ]
 [-2.067       3.58014902  0.        ]
 [ 0.          0.          6.652     ]]
Number of wannier functions: 13
Number of R points: 267
Recommended size of FFT grid [6 6 4]
setting Rvec
setting AA..
setting AA - OK
Real-space lattice:
 [[ 4.134       0.          0.        ]
 [-2.067       3.58014902  0.        ]
 [ 0.          0.          6.652     ]]
Number of wannier functions: 13
Number of R points: 267
Recommended size of FFT grid [6 6 4]
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/utility.py:68: UserWarning: usually need to provide either with real or reciprocal lattice.If you only want to generate a random symmetric tensor - that it fine
  warnings.warn("usually need to provide either with real or reciprocal lattice."
nspin in SOC: 2
setting the Rvectors with wannier centers (cart):
 [[ 0.00000000e+00  0.00000000e+00  0.00000000e+00]
 [ 0.00000000e+00  0.00000000e+00  3.32600000e+00]
 [ 0.00000000e+00  0.00000000e+00  0.00000000e+00]
 [ 0.00000000e+00  0.00000000e+00  3.32600000e+00]
 [ 0.00000000e+00  0.00000000e+00  0.00000000e+00]
 [ 0.00000000e+00  0.00000000e+00  3.32600000e+00]
 [ 0.00000000e+00  0.00000000e+00  0.00000000e+00]
 [ 0.00000000e+00  0.00000000e+00  3.32600000e+00]
 [ 0.00000000e+00  0.00000000e+00  0.00000000e+00]
 [ 0.00000000e+00  0.00000000e+00  3.32600000e+00]
 [ 7.05488710e-01  2.79408011e+00  1.63112851e+00]
 [ 7.05546670e-01  2.79411357e+00  1.69453500e+00]
 [-7.05488710e-01  2.79408011e+00  1.63112851e+00]
 [-7.05546670e-01  2.79411357e+00  1.69453500e+00]
 [ 2.00006889e-16  1.57213782e+00  1.63112851e+00]
 [-1.16724460e-17  1.57207089e+00  1.69453500e+00]
 [ 1.36151129e+00  7.86068910e-01  5.02087149e+00]
 [ 1.36145333e+00  7.86035446e-01  4.95746500e+00]
 [ 2.77248871e+00  7.86068910e-01  5.02087149e+00]
 [ 2.77254667e+00  7.86035446e-01  4.95746500e+00]
 [ 2.06700000e+00  2.00801120e+00  5.02087149e+00]
 [ 2.06700000e+00  2.00807813e+00  4.95746500e+00]
 [ 1.98359847e-17  2.38676601e+00  1.92689851e+00]
 [ 1.98359847e-17  2.38676601e+00  1.40012024e+00]
 [ 2.06700000e+00  1.19338301e+00  4.72510149e+00]
 [ 2.06700000e+00  1.19338301e+00  5.25187976e+00]]
 (red):
 [[0.         0.         0.        ]
 [0.         0.         0.5       ]
 [0.         0.         0.        ]
 [0.         0.         0.5       ]
 [0.         0.         0.        ]
 [0.         0.         0.5       ]
 [0.         0.         0.        ]
 [0.         0.         0.5       ]
 [0.         0.         0.        ]
 [0.         0.         0.5       ]
 [0.56087364 0.78043682 0.24520874]
 [0.56089233 0.78044617 0.25474068]
 [0.21956318 0.78043682 0.24520874]
 [0.21955383 0.78044617 0.25474068]
 [0.21956318 0.43912636 0.24520874]
 [0.21955383 0.43910767 0.25474068]
 [0.43912636 0.21956318 0.75479126]
 [0.43910767 0.21955383 0.74525932]
 [0.78043682 0.21956318 0.75479126]
 [0.78044617 0.21955383 0.74525932]
 [0.78043682 0.56087364 0.75479126]
 [0.78044617 0.56089233 0.74525932]
 [0.33333333 0.66666667 0.28967206]
 [0.33333333 0.66666667 0.2104811 ]
 [0.66666667 0.33333333 0.71032794]
 [0.66666667 0.33333333 0.7895189 ]]
setting Rvec
using magmoms
 [[ 0.00000000e+00  0.00000000e+00  4.40436811e+00]
 [-0.00000000e+00 -0.00000000e+00 -4.40436811e+00]
 [-0.00000000e+00 -0.00000000e+00 -2.92436898e-10]
 [-0.00000000e+00 -0.00000000e+00 -2.92436898e-10]]
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sym_wann_2.py:258: UserWarning: symmetrization of matrices {'dV_soc_wann_0_0'} is not tested. use on your own risk
  warnings.warn(f"symmetrization of matrices {unknown} is not tested. use on your own risk")
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sym_wann_2.py:258: UserWarning: symmetrization of matrices {'dV_soc_wann_1_1'} is not tested. use on your own risk
  warnings.warn(f"symmetrization of matrices {unknown} is not tested. use on your own risk")
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/symmetry/sym_wann_2.py:258: UserWarning: symmetrization of matrices {'overlap_up_down', 'dV_soc_wann_0_1'} is not tested. use on your own risk
  warnings.warn(f"symmetrization of matrices {unknown} is not tested. use on your own risk")
using magmoms
 [[ 0.00000000e+00  0.00000000e+00  4.40436811e+00]
 [-0.00000000e+00 -0.00000000e+00 -4.40436811e+00]
 [-0.00000000e+00 -0.00000000e+00 -2.92436898e-10]
 [-0.00000000e+00 -0.00000000e+00 -2.92436898e-10]]
using magmoms
 [[ 2.69689765e-16  4.40436811e+00  2.69689765e-16]
 [-2.69689765e-16 -4.40436811e+00 -2.69689765e-16]
 [-1.79065955e-26 -2.92436898e-10 -1.79065955e-26]
 [-1.79065955e-26 -2.92436898e-10 -1.79065955e-26]]
Saving SystemSOC to system_soc
saving system of class SystemSOC to system_soc
 properties: ['num_wann', 'real_lattice', 'iRvec', 'periodic', 'is_phonon', 'wannier_centers_cart', 'pointgroup', 'cell']
saving num_wann to system_soc/num_wann.npz
saving real_lattice to system_soc/real_lattice.npz
saving iRvec to system_soc/iRvec.npz
saving periodic to system_soc/periodic.npz
saving is_phonon to system_soc/is_phonon.npz
saving wannier_centers_cart to system_soc/wannier_centers_cart.npz
saving pointgroup to system_soc/pointgroup.npz
saving cell to system_soc/cell.npz
saving system of class System_w90 to system_soc/system_up
 properties: ['num_wann', 'real_lattice', 'iRvec', 'periodic', 'is_phonon', 'wannier_centers_cart', 'pointgroup']
saving num_wann
saving num_wann to system_soc/system_up/num_wann.npz
 - Ok!
saving real_lattice
saving real_lattice to system_soc/system_up/real_lattice.npz
 - Ok!
saving iRvec
saving iRvec to system_soc/system_up/iRvec.npz
 - Ok!
saving periodic
saving periodic to system_soc/system_up/periodic.npz
 - Ok!
saving is_phonon
saving is_phonon to system_soc/system_up/is_phonon.npz
 - Ok!
saving wannier_centers_cart
saving wannier_centers_cart to system_soc/system_up/wannier_centers_cart.npz
 - Ok!
saving pointgroup
saving pointgroup to system_soc/system_up/pointgroup.npz
 - Ok!
saving Ham - Ok!
saving AA - Ok!
saving system of class System_w90 to system_soc/system_down
 properties: ['num_wann', 'real_lattice', 'iRvec', 'periodic', 'is_phonon', 'wannier_centers_cart', 'pointgroup']
saving num_wann
saving num_wann to system_soc/system_down/num_wann.npz
 - Ok!
saving real_lattice
saving real_lattice to system_soc/system_down/real_lattice.npz
 - Ok!
saving iRvec
saving iRvec to system_soc/system_down/iRvec.npz
 - Ok!
saving periodic
saving periodic to system_soc/system_down/periodic.npz
 - Ok!
saving is_phonon
saving is_phonon to system_soc/system_down/is_phonon.npz
 - Ok!
saving wannier_centers_cart
saving wannier_centers_cart to system_soc/system_down/wannier_centers_cart.npz
 - Ok!
saving pointgroup
saving pointgroup to system_soc/system_down/pointgroup.npz
 - Ok!
saving Ham - Ok!
saving AA - Ok!

Compute the wannierized bandstructure with SOC along the high-symmetry path

[10]:
from wannierberri.grid import Path
from wannierberri.evaluate_k import evaluate_k_path
from wannierberri.system.system_soc import SystemSOC
system_soc = SystemSOC.from_npz("system_soc")
kz = 0.35 / (2*np.pi/c)
path = Path(lattice,
            nodes=[
                [2 / 3, -1 / 3, 0],
                [0, 0, 0],
                [-2 / 3, 1 / 3, 0],
                None,
                [-0.5, 0, kz],
                [0, 0, kz],
                [0.5, 0, kz],
            ],
            labels=[r"${\rm K}\leftarrow$",
                    r"$\Gamma$",
                    r"$\rightarrow{\rm K}$",
                    r"$\overline{\rm M}\leftarrow$",
                    r"$\overline{\Gamma}$",
                    r"$\rightarrow\overline{\rm M}$"],
            length=300)   # length [ Ang] ~= 2*pi/dk

bands_wannier_soc= evaluate_k_path(system_soc, path=path, quantities=["spin"], parallel=parallel_env)
bands_wannier_up = evaluate_k_path(system_soc.system_up, path=path, parallel=parallel_env)
bands_wannier_dw = evaluate_k_path(system_soc.system_down, path=path, parallel=parallel_env)
loading real_lattice - Ok!
loading wannier_centers_cart - Ok!
loading pointgroup - Ok!
loading iRvec - Ok!
loading periodic - Ok!
loading is_phonon - Ok!
loading num_wann - Ok!
loading R_matrix Ham - Ok!
loading R_matrix AA - Ok!
loading real_lattice - Ok!
loading wannier_centers_cart - Ok!
loading pointgroup - Ok!
loading iRvec - Ok!
loading periodic - Ok!
loading is_phonon - Ok!
loading num_wann - Ok!
loading R_matrix Ham - Ok!
loading R_matrix AA - Ok!
Starting run()
Using the follwing calculators :
############################################################

 'tabulate'  :  <wannierberri.calculators.tabulate.TabulatorAll object at 0x74c0905060c0>  :
    TabulatorAll - a pack of all k-resolved calculators (Tabulators)

 Includes the following tabulators :
--------------------------------------------------
 "spin" : <wannierberri.calculators.tabulate.Spin object at 0x74c0e9e8a540> :  Spin expectation :math:` \langle u | \mathbf{\sigma} | u \rangle`

 "Energy" : <wannierberri.calculators.tabulate.Energy object at 0x74c0606b6e70> : calculator not described

--------------------------------------------------

############################################################
Calculation along a path - checking calculators for compatibility
tabulate <wannierberri.calculators.tabulate.TabulatorAll object at 0x74c0905060c0>
All calculators are compatible
Symmetrization switched off for Path
Grid is regular
The set of k points is a Path() with 181 points and labels {0: '${\\rm K}\\leftarrow$', 48: '$\\Gamma$', 96: '$\\overline{\\rm M}\\leftarrow$', 138: '$\\overline{\\Gamma}$', 180: '$\\rightarrow\\overline{\\rm M}$'}
generating K_list
Done
Done, sum of weights:181.0

############################################################
Iteration 0 out of 0
processing 181 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/utility.py:68: UserWarning: usually need to provide either with real or reciprocal lattice.If you only want to generate a random symmetric tensor - that it fine
  warnings.warn("usually need to provide either with real or reciprocal lattice."
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/grid/path.py:163: UserWarning: symmetry is not used for a tabulation along path
  warnings.warn("symmetry is not used for a tabulation along path")
time for processing    181 K-points on  16 processes:     0.3538 ; per K-point          0.0020 ; proc-sec per K-point          0.0313
time1 =  0.0051348209381103516
Totally processed 181 K-points
run() finished
Starting run()
Using the follwing calculators :
############################################################

 'tabulate'  :  <wannierberri.calculators.tabulate.TabulatorAll object at 0x74c0905060c0>  :
    TabulatorAll - a pack of all k-resolved calculators (Tabulators)

 Includes the following tabulators :
--------------------------------------------------
 "Energy" : <wannierberri.calculators.tabulate.Energy object at 0x74c0906187d0> : calculator not described

--------------------------------------------------

############################################################
Calculation along a path - checking calculators for compatibility
tabulate <wannierberri.calculators.tabulate.TabulatorAll object at 0x74c0905060c0>
All calculators are compatible
Symmetrization switched off for Path
Grid is regular
The set of k points is a Path() with 181 points and labels {0: '${\\rm K}\\leftarrow$', 48: '$\\Gamma$', 96: '$\\overline{\\rm M}\\leftarrow$', 138: '$\\overline{\\Gamma}$', 180: '$\\rightarrow\\overline{\\rm M}$'}
generating K_list
Done
Done, sum of weights:181.0

############################################################
Iteration 0 out of 0
processing 181 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing    181 K-points on  16 processes:     0.0896 ; per K-point          0.0005 ; proc-sec per K-point          0.0079
time1 =  0.003645658493041992
Totally processed 181 K-points
run() finished
Starting run()
Using the follwing calculators :
############################################################

 'tabulate'  :  <wannierberri.calculators.tabulate.TabulatorAll object at 0x74c0905060c0>  :
    TabulatorAll - a pack of all k-resolved calculators (Tabulators)

 Includes the following tabulators :
--------------------------------------------------
 "Energy" : <wannierberri.calculators.tabulate.Energy object at 0x74c090662e40> : calculator not described

--------------------------------------------------

############################################################
Calculation along a path - checking calculators for compatibility
tabulate <wannierberri.calculators.tabulate.TabulatorAll object at 0x74c0905060c0>
All calculators are compatible
Symmetrization switched off for Path
Grid is regular
The set of k points is a Path() with 181 points and labels {0: '${\\rm K}\\leftarrow$', 48: '$\\Gamma$', 96: '$\\overline{\\rm M}\\leftarrow$', 138: '$\\overline{\\Gamma}$', 180: '$\\rightarrow\\overline{\\rm M}$'}
generating K_list
Done
Done, sum of weights:181.0

############################################################
Iteration 0 out of 0
processing 181 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing    181 K-points on  16 processes:     0.0818 ; per K-point          0.0005 ; proc-sec per K-point          0.0072
time1 =  0.0038216114044189453
Totally processed 181 K-points
run() finished

Plot the Bandstructure with and without SOC

[12]:
from matplotlib import pyplot as plt

EF = 6.885145845031927

fig, axes = plt.subplots(4, 1, sharey=True, sharex=True, figsize=(15, 15))

for i in range(3):
    component = "xyz"[i]
    bands_wannier_soc.plot_path_fat(path=path,
                       Eshift=EF,
                       quantity="spin",
                       component=component,
                       mode="color",
                       label=f"soc_nscf, S{component}, th={theta}, phi={phi}",
                        axes=axes[1 + i],
                       fatmax=10,
                        linecolor="orange",
                        close_fig=False,
                        show_fig=False,
                        kwargs_line=dict(linestyle='-', lw=0.0),
    )


bands_wannier_up.plot_path_fat(path=path,
                       Eshift=EF,
                       axes=axes[0],
                       label="spin-up",
                       linecolor="red",
                       close_fig=False,
                       show_fig=False,
)

bands_wannier_dw.plot_path_fat(path=path,
                       label="spin-dw",
                       axes=axes[0],
                       Eshift=EF,
                       linecolor="blue",
                       close_fig=False,
                       show_fig=False,)

plt.ylim(-2,0)
plt.show()

../../../_images/tutorials_8_GPAW_3.MnTe_mnte_22_0.png

AHC

magnetic moments along “y” direction

Below we evaluate only the “internal” (which depend only on the Hamiltonian) contributions to the anomalous Hall conductivity (AHC), how to include the “external” contributions is explained in the previous tutorial on bcc Fe.

[13]:
import os
import numpy as np
import wannierberri as wb
from wannierberri.system.system_soc import SystemSOC


system_soc_y = SystemSOC.from_npz("system_soc")
theta=90
phi=90
system_soc_y.set_soc_axis(theta=theta, phi=phi, units="degrees")

grid = wb.grid.Grid(system_soc_y, NK=200)

EF = 6.885145845031927
tetra = False

Efermi = np.linspace(EF - 2, EF, 1001)

os.makedirs("results", exist_ok=True)
result_int = wb.run(system_soc_y,
        grid=grid,
        parallel=parallel_env,
        fout_name=f"results/{seed}-soc-y",
        calculators={"ahc_int": wb.calculators.static.AHC(Efermi=Efermi, tetra=tetra, kwargs_formula={"external_terms": False})},
        adpt_num_iter=50,
        restart=False,
        print_progress_step=5,
        )

loading real_lattice - Ok!
loading wannier_centers_cart - Ok!
loading pointgroup - Ok!
loading iRvec - Ok!
loading periodic - Ok!
loading is_phonon - Ok!
loading num_wann - Ok!
loading R_matrix Ham - Ok!
loading R_matrix AA - Ok!
loading real_lattice - Ok!
loading wannier_centers_cart - Ok!
loading pointgroup - Ok!
loading iRvec - Ok!
loading periodic - Ok!
loading is_phonon - Ok!
loading num_wann - Ok!
loading R_matrix Ham - Ok!
loading R_matrix AA - Ok!
/home/stepan/github/wannier-berri-work/wannier-berri-soc-nscf/wannierberri/utility.py:68: UserWarning: usually need to provide either with real or reciprocal lattice.If you only want to generate a random symmetric tensor - that it fine
  warnings.warn("usually need to provide either with real or reciprocal lattice."
using magmoms
 [[ 2.69689765e-16  4.40436811e+00  2.69689765e-16]
 [-2.69689765e-16 -4.40436811e+00 -2.69689765e-16]
 [-1.79065955e-26 -2.92436898e-10 -1.79065955e-26]
 [-1.79065955e-26 -2.92436898e-10 -1.79065955e-26]]
Minimal symmetric FFT grid :  [11 11  7]
Starting run()
Using the follwing calculators :
############################################################

 'ahc_int'  :  <wannierberri.calculators.static.AHC object at 0x74c0e06bb140>  : Anomalous Hall conductivity (:math:`s^3 \cdot A^2 / (kg \cdot m^3) = S/m`)

        | With Fermi sea integral Eq(11) in `Ref <https://www.nature.com/articles/s41524-021-00498-5>`__
        | Output: :math:`O = -e^2/\hbar \int [dk] \Omega f`
        | Instruction: :math:`j_\alpha = \sigma_{\alpha\beta} E_\beta = \epsilon_{\alpha\beta\delta} O_\delta E_\beta`
############################################################
Calculation on  grid - checking calculators for compatibility
ahc_int <wannierberri.calculators.static.AHC object at 0x74c0e06bb140>
All calculators are compatible
Grid is regular
The set of k points is a Grid() with NKdiv=[10 10 25], NKFFT=[20 20  8], NKtot=[200 200 200]
generating K_list
Done in 0.012741565704345703 s
excluding symmetry-equivalent K-points from initial grid
Done in 0.04427790641784668 s
K_list contains 403 Irreducible points(16.12%) out of initial 10x10x25=2500 grid
Done, sum of weights:0.999999999999999

############################################################
Iteration 0 out of 50
processing 403 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.2                 173.8                 180.9
                  48             15.0                 111.2                 126.2
                  80             24.1                  97.5                 121.6
                 112             32.8                  85.3                 118.1
                 144             41.4                  74.4                 115.8
                 176             49.7                  64.0                 113.7
                 208             58.1                  54.5                 112.7
                 240             66.9                  45.4                 112.3
                 272             75.6                  36.4                 112.0
                 304             84.0                  27.4                 111.4
                 336             92.9                  18.5                 111.4
                 368            101.4                   9.6                 111.1
                 400            109.4                   0.8                 110.2
time for processing    403 K-points on  16 processes:   109.8644 ; per K-point          0.2726 ; proc-sec per K-point          4.3619
time1 =  0.008203744888305664
time2 =  0.011982202529907227
sum of weights now :0.999999999999999

############################################################
Iteration 1 out of 50
processing 4 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      4 K-points on  16 processes:     2.6180 ; per K-point          0.6545 ; proc-sec per K-point         10.4722
time1 =  0.00013756752014160156
time2 =  0.006268501281738281
sum of weights now :0.9999999999999989

############################################################
Iteration 2 out of 50
processing 4 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      4 K-points on  16 processes:     2.5903 ; per K-point          0.6476 ; proc-sec per K-point         10.3612
time1 =  0.0002353191375732422
time2 =  0.006016254425048828
sum of weights now :0.9999999999999987

############################################################
Iteration 3 out of 50
processing 8 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      8 K-points on  16 processes:     3.6854 ; per K-point          0.4607 ; proc-sec per K-point          7.3707
time1 =  0.00021409988403320312
time2 =  0.006436824798583984
sum of weights now :0.9999999999999987

############################################################
Iteration 4 out of 50
processing 4 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      4 K-points on  16 processes:     2.6711 ; per K-point          0.6678 ; proc-sec per K-point         10.6842
time1 =  0.0001742839813232422
time2 =  0.006453752517700195
sum of weights now :0.9999999999999986

############################################################
Iteration 5 out of 50
processing 4 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      4 K-points on  16 processes:     2.4979 ; per K-point          0.6245 ; proc-sec per K-point          9.9918
time1 =  0.0002307891845703125
time2 =  0.007687807083129883
sum of weights now :0.9999999999999984

############################################################
Iteration 6 out of 50
processing 4 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      4 K-points on  16 processes:     2.5884 ; per K-point          0.6471 ; proc-sec per K-point         10.3535
time1 =  0.0002791881561279297
time2 =  0.007031917572021484
sum of weights now :0.9999999999999983

############################################################
Iteration 7 out of 50
processing 8 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      8 K-points on  16 processes:     3.4217 ; per K-point          0.4277 ; proc-sec per K-point          6.8434
time1 =  0.00016570091247558594
time2 =  0.005188465118408203
sum of weights now :0.9999999999999982

############################################################
Iteration 8 out of 50
processing 4 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      4 K-points on  16 processes:     2.4196 ; per K-point          0.6049 ; proc-sec per K-point          9.6786
time1 =  0.00021219253540039062
time2 =  0.005506992340087891
sum of weights now :0.9999999999999981

############################################################
Iteration 9 out of 50
processing 4 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      4 K-points on  16 processes:     2.5337 ; per K-point          0.6334 ; proc-sec per K-point         10.1348
time1 =  0.0002913475036621094
time2 =  0.005382537841796875
sum of weights now :0.999999999999998

############################################################
Iteration 10 out of 50
processing 4 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      4 K-points on  16 processes:     2.3921 ; per K-point          0.5980 ; proc-sec per K-point          9.5685
time1 =  0.0002586841583251953
time2 =  0.006248950958251953
sum of weights now :0.999999999999998

############################################################
Iteration 11 out of 50
processing 4 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      4 K-points on  16 processes:     2.3319 ; per K-point          0.5830 ; proc-sec per K-point          9.3277
time1 =  0.00035119056701660156
time2 =  0.006663322448730469
sum of weights now :0.9999999999999974

############################################################
Iteration 12 out of 50
processing 12 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     12 K-points on  16 processes:     4.8746 ; per K-point          0.4062 ; proc-sec per K-point          6.4995
time1 =  0.0006597042083740234
time2 =  0.006142854690551758
sum of weights now :0.999999999999997

############################################################
Iteration 13 out of 50
processing 12 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     12 K-points on  16 processes:     4.7689 ; per K-point          0.3974 ; proc-sec per K-point          6.3586
time1 =  0.0006327629089355469
time2 =  0.01129603385925293
sum of weights now :0.999999999999997

############################################################
Iteration 14 out of 50
processing 8 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      8 K-points on  16 processes:     3.4274 ; per K-point          0.4284 ; proc-sec per K-point          6.8548
time1 =  0.0001983642578125
time2 =  0.0052947998046875
sum of weights now :0.9999999999999968

############################################################
Iteration 15 out of 50
processing 8 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing      8 K-points on  16 processes:     3.3783 ; per K-point          0.4223 ; proc-sec per K-point          6.7566
time1 =  0.00046062469482421875
time2 =  0.009725809097290039
sum of weights now :0.9999999999999963

############################################################
Iteration 16 out of 50
processing 16 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     16 K-points on  16 processes:     5.8587 ; per K-point          0.3662 ; proc-sec per K-point          5.8587
time1 =  0.0003521442413330078
time2 =  0.006136655807495117
sum of weights now :0.9999999999999958

############################################################
Iteration 17 out of 50
processing 12 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     12 K-points on  16 processes:     4.9260 ; per K-point          0.4105 ; proc-sec per K-point          6.5681
time1 =  0.0002868175506591797
time2 =  0.0057871341705322266
sum of weights now :0.999999999999995

############################################################
Iteration 18 out of 50
processing 20 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              6.5                   1.6                   8.2
time for processing     20 K-points on  16 processes:     6.9187 ; per K-point          0.3459 ; proc-sec per K-point          5.5350
time1 =  0.00044608116149902344
time2 =  0.007273197174072266
sum of weights now :0.9999999999999941

############################################################
Iteration 19 out of 50
processing 20 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              6.8                   1.7                   8.5
time for processing     20 K-points on  16 processes:     7.0971 ; per K-point          0.3549 ; proc-sec per K-point          5.6777
time1 =  0.0005843639373779297
time2 =  0.006464242935180664
sum of weights now :0.9999999999999936

############################################################
Iteration 20 out of 50
processing 18 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              6.2                   0.8                   7.0
time for processing     18 K-points on  16 processes:     6.6277 ; per K-point          0.3682 ; proc-sec per K-point          5.8913
time1 =  0.0005707740783691406
time2 =  0.01274728775024414
sum of weights now :0.9999999999999931

############################################################
Iteration 21 out of 50
processing 12 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     12 K-points on  16 processes:     4.8030 ; per K-point          0.4002 ; proc-sec per K-point          6.4040
time1 =  0.0003199577331542969
time2 =  0.006445407867431641
sum of weights now :0.9999999999999926

############################################################
Iteration 22 out of 50
processing 16 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     16 K-points on  16 processes:     5.7218 ; per K-point          0.3576 ; proc-sec per K-point          5.7218
time1 =  0.00040841102600097656
time2 =  0.00537419319152832
sum of weights now :0.9999999999999921

############################################################
Iteration 23 out of 50
processing 12 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     12 K-points on  16 processes:     4.9905 ; per K-point          0.4159 ; proc-sec per K-point          6.6541
time1 =  0.0006797313690185547
time2 =  0.008481502532958984
sum of weights now :0.9999999999999916

############################################################
Iteration 24 out of 50
processing 20 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              6.6                   1.7                   8.3
time for processing     20 K-points on  16 processes:     7.1416 ; per K-point          0.3571 ; proc-sec per K-point          5.7132
time1 =  0.0005571842193603516
time2 =  0.0065593719482421875
sum of weights now :0.9999999999999908

############################################################
Iteration 25 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.0                   3.5                  10.5
time for processing     24 K-points on  16 processes:     7.9874 ; per K-point          0.3328 ; proc-sec per K-point          5.3250
time1 =  0.0005202293395996094
time2 =  0.006932735443115234
sum of weights now :0.9999999999999898

############################################################
Iteration 26 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.0                   3.5                  10.5
time for processing     24 K-points on  16 processes:     7.9818 ; per K-point          0.3326 ; proc-sec per K-point          5.3212
time1 =  0.0005197525024414062
time2 =  0.006588459014892578
sum of weights now :0.9999999999999893

############################################################
Iteration 27 out of 50
processing 20 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              6.3                   1.6                   7.9
time for processing     20 K-points on  16 processes:     6.7701 ; per K-point          0.3385 ; proc-sec per K-point          5.4160
time1 =  0.00036597251892089844
time2 =  0.005646705627441406
sum of weights now :0.9999999999999889

############################################################
Iteration 28 out of 50
processing 16 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     16 K-points on  16 processes:     5.5389 ; per K-point          0.3462 ; proc-sec per K-point          5.5389
time1 =  0.0004162788391113281
time2 =  0.006862640380859375
sum of weights now :0.9999999999999883

############################################################
Iteration 29 out of 50
processing 20 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              6.7                   1.7                   8.3
time for processing     20 K-points on  16 processes:     7.2080 ; per K-point          0.3604 ; proc-sec per K-point          5.7664
time1 =  0.0007338523864746094
time2 =  0.00673222541809082
sum of weights now :0.9999999999999877

############################################################
Iteration 30 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.0                   3.5                  10.5
time for processing     24 K-points on  16 processes:     8.2766 ; per K-point          0.3449 ; proc-sec per K-point          5.5177
time1 =  0.0005040168762207031
time2 =  0.006816387176513672
sum of weights now :0.9999999999999872

############################################################
Iteration 31 out of 50
processing 20 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              6.8                   1.7                   8.5
time for processing     20 K-points on  16 processes:     7.2036 ; per K-point          0.3602 ; proc-sec per K-point          5.7629
time1 =  0.00047016143798828125
time2 =  0.00696563720703125
sum of weights now :0.999999999999987

############################################################
Iteration 32 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.3                   3.6                  10.9
time for processing     24 K-points on  16 processes:     8.4053 ; per K-point          0.3502 ; proc-sec per K-point          5.6035
time1 =  0.000659942626953125
time2 =  0.008317232131958008
sum of weights now :0.999999999999986

############################################################
Iteration 33 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.3                   3.6                  10.9
time for processing     24 K-points on  16 processes:     7.9262 ; per K-point          0.3303 ; proc-sec per K-point          5.2841
time1 =  0.0004787445068359375
time2 =  0.006205320358276367
sum of weights now :0.9999999999999852

############################################################
Iteration 34 out of 50
processing 20 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              6.8                   1.7                   8.5
time for processing     20 K-points on  16 processes:     7.2445 ; per K-point          0.3622 ; proc-sec per K-point          5.7956
time1 =  0.0011258125305175781
time2 =  0.008062362670898438
sum of weights now :0.9999999999999848

############################################################
Iteration 35 out of 50
processing 16 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     16 K-points on  16 processes:     6.3314 ; per K-point          0.3957 ; proc-sec per K-point          6.3314
time1 =  0.00034999847412109375
time2 =  0.006597757339477539
sum of weights now :0.9999999999999841

############################################################
Iteration 36 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.4                   3.7                  11.0
time for processing     24 K-points on  16 processes:     8.3841 ; per K-point          0.3493 ; proc-sec per K-point          5.5894
time1 =  0.0008716583251953125
time2 =  0.008568763732910156
sum of weights now :0.9999999999999832

############################################################
Iteration 37 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.3                   3.7                  11.0
time for processing     24 K-points on  16 processes:     8.3491 ; per K-point          0.3479 ; proc-sec per K-point          5.5661
time1 =  0.0004677772521972656
time2 =  0.0069620609283447266
sum of weights now :0.9999999999999823

############################################################
Iteration 38 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.4                   3.7                  11.1
time for processing     24 K-points on  16 processes:     8.2163 ; per K-point          0.3423 ; proc-sec per K-point          5.4775
time1 =  0.0005958080291748047
time2 =  0.008387565612792969
sum of weights now :0.9999999999999816

############################################################
Iteration 39 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.4                   3.7                  11.1
time for processing     24 K-points on  16 processes:     8.1330 ; per K-point          0.3389 ; proc-sec per K-point          5.4220
time1 =  0.0010483264923095703
time2 =  0.01163172721862793
sum of weights now :0.999999999999981

############################################################
Iteration 40 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.2                   3.6                  10.8
time for processing     24 K-points on  16 processes:     8.0699 ; per K-point          0.3362 ; proc-sec per K-point          5.3799
time1 =  0.0006461143493652344
time2 =  0.007129669189453125
sum of weights now :0.9999999999999807

############################################################
Iteration 41 out of 50
processing 16 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     16 K-points on  16 processes:     6.0274 ; per K-point          0.3767 ; proc-sec per K-point          6.0274
time1 =  0.0004124641418457031
time2 =  0.008254289627075195
sum of weights now :0.9999999999999803

############################################################
Iteration 42 out of 50
processing 16 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     16 K-points on  16 processes:     6.0304 ; per K-point          0.3769 ; proc-sec per K-point          6.0304
time1 =  0.0007953643798828125
time2 =  0.006538867950439453
sum of weights now :0.99999999999998

############################################################
Iteration 43 out of 50
processing 20 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              6.9                   1.7                   8.6
time for processing     20 K-points on  16 processes:     7.3302 ; per K-point          0.3665 ; proc-sec per K-point          5.8642
time1 =  0.00044608116149902344
time2 =  0.007014751434326172
sum of weights now :0.9999999999999793

############################################################
Iteration 44 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.3                   3.7                  11.0
time for processing     24 K-points on  16 processes:     7.9420 ; per K-point          0.3309 ; proc-sec per K-point          5.2947
time1 =  0.00047850608825683594
time2 =  0.006898641586303711
sum of weights now :0.9999999999999787

############################################################
Iteration 45 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.0                   3.5                  10.5
time for processing     24 K-points on  16 processes:     8.1079 ; per K-point          0.3378 ; proc-sec per K-point          5.4052
time1 =  0.0013275146484375
time2 =  0.010807991027832031
sum of weights now :0.9999999999999777

############################################################
Iteration 46 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.2                   3.6                  10.9
time for processing     24 K-points on  16 processes:     8.2289 ; per K-point          0.3429 ; proc-sec per K-point          5.4859
time1 =  0.0006849765777587891
time2 =  0.007607221603393555
sum of weights now :0.9999999999999768

############################################################
Iteration 47 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.3                   3.7                  11.0
time for processing     24 K-points on  16 processes:     7.9600 ; per K-point          0.3317 ; proc-sec per K-point          5.3067
time1 =  0.0005161762237548828
time2 =  0.007196187973022461
sum of weights now :0.9999999999999765

############################################################
Iteration 48 out of 50
processing 16 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     16 K-points on  16 processes:     5.9199 ; per K-point          0.3700 ; proc-sec per K-point          5.9199
time1 =  0.00036978721618652344
time2 =  0.006857872009277344
sum of weights now :0.9999999999999761

############################################################
Iteration 49 out of 50
processing 16 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
time for processing     16 K-points on  16 processes:     5.9352 ; per K-point          0.3710 ; proc-sec per K-point          5.9352
time1 =  0.00040435791015625
time2 =  0.006731748580932617
sum of weights now :0.9999999999999751

############################################################
Iteration 50 out of 50
processing 24 K points : using  16 processes.
# K-points calculated  Wall time (sec)  Est. remaining (sec)   Est. total (sec)
                  16              7.2                   3.6                  10.8
time for processing     24 K-points on  16 processes:     8.0830 ; per K-point          0.3368 ; proc-sec per K-point          5.3887
time1 =  0.0014445781707763672
Totally processed 1197 K-points
run() finished

plot AHC vs energy

[34]:
from matplotlib import pyplot as plt
import numpy as np
theta=0
phi=0


fig, axes = plt.subplots(1, 1, sharey=False, sharex=True, figsize=(10, 8))
axes=[axes]
for i in range(0,51,5):
    color = (0, 1-i/50, i/50)
    res =  np.load(f"results/MnTe-soc-y-ahc_int_iter-{i:04d}.npz")
    Efermi = res['Energies_0']
    ahc = res["data"]
    axes[0].plot(Efermi - EF, ahc[:,2]/100, label=f"iter={i}", color=color)
axes[0].set_title("AHC internal terms convergence")
# axes[0].set_ylim(-1000, 0)
axes[0].set_xlabel("E-EF (eV)")
axes[0].set_ylabel("AHC (S/cm)")
axes[0].legend()
plt.tight_layout()
plt.savefig("ahc_convergence.png", dpi=300)
../../../_images/tutorials_8_GPAW_3.MnTe_mnte_26_0.png

We can observe 3 things, which are generally true, not only for AHC, but for most other effects:

  1. The internal terms converge slower than the external ones.

  2. The external terms usually give a smaller contribution. Often (but not always) it is negligible. However, check it for a particular material and effect.

  3. The external terms take more time to evaluate per k-point.

Thus, it makes sense to evaluate them separately, using different k-point grids and number of iterations for the adaptive integration.

Other directions of magnetic moments

We could try to calculate the AHC when the magnetic moments are along “z” or “x” directions. However, we can see that in those cases the point group contains symmetries which force the AHC to be zero. This can be demonstrated by the following way:

[38]:
from irrep.spacegroup import SpaceGroup
from wannierberri.symmetry.point_symmetry import PointGroup

#  x direction
mg = SpaceGroup.from_cell(real_lattice=lattice, positions=positions, typat=[1,1,2,2], magmom=[[1,0,0],[-1,0,0],[0,0,0],[0,0,0]])
# mg.show()
pointgroup = PointGroup(spacegroup=mg)
print (f"components of AHC for magnetic moments along x direction: {pointgroup.get_symmetric_components(rank=1, TRodd=True, Iodd=False)}")

mg = SpaceGroup.from_cell(real_lattice=lattice, positions=positions, typat=[1,1,2,2], magmom=[[0,1,0],[0,-1,0],[0,0,0],[0,0,0]])
# mg.show()
pointgroup = PointGroup(spacegroup=mg)
print (f"components of AHC for magnetic moments along y direction: {pointgroup.get_symmetric_components(rank=1, TRodd=True, Iodd=False)}")

mg = SpaceGroup.from_cell(real_lattice=lattice, positions=positions, typat=[1,1,2,2], magmom=[[0,0,1],[0,0,-1],[0,0,0],[0,0,0]])
# mg.show()
pointgroup = PointGroup(spacegroup=mg)
print (f"components of AHC for magnetic moments along z direction: {pointgroup.get_symmetric_components(rank=1, TRodd=True, Iodd=False)}")
components of AHC for magnetic moments along x direction: ['0=x=y=z']
components of AHC for magnetic moments along y direction: ['0=x=y', 'z']
components of AHC for magnetic moments along z direction: ['0=x=y=z']