Skip to content

API Reference

wannier90io.parse_win_raw(string: str) -> dict

Parse WIN

Parameters:

Name Type Description Default
string str

the WIN text

required

Returns:

Type Description
dict

the parsed WIN

Source code in wannier90io/_win.py
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
def parse_win_raw(string: str) -> dict:
    """
    Parse WIN

    Arguments:
        string: the WIN text

    Returns:
        the parsed WIN
    """
    comments = _core.extract_comments(string)
    parameters = _core.parse_parameters(_core.extract_parameters(string))
    blocks = _core.parse_blocks(_core.extract_blocks(string))

    parsed_win = {
        'comments': comments,
        'parameters': parameters,
        'blocks': blocks,
    }
    if 'unit_cell_cart' in blocks:
        parsed_win['unit_cell_cart'] = parse_unit_cell(blocks['unit_cell_cart'])
    if 'atoms_cart' in blocks:
        parsed_win['atoms_cart'] = parse_atoms(blocks['atoms_cart'])
    if 'atoms_frac' in blocks:
        parsed_win['atoms_frac'] = parse_atoms(blocks['atoms_frac'])
    if 'projections' in blocks:
        parsed_win['projections'] = parse_projections(blocks['projections'])
    if 'kpoints' in blocks:
        parsed_win['kpoints'] = parse_kpoints(blocks['kpoints'])

    return parsed_win

wannier90io.parse_nnkp_raw(string: str) -> dict

Parse NNKP

Parameters:

Name Type Description Default
string str

the NNKP text

required

Returns:

Type Description
dict

the parsed NNKP

Source code in wannier90io/_nnkp.py
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
def parse_nnkp_raw(string: str) -> dict:
    """
    Parse NNKP

    Arguments:
        string: the NNKP text

    Returns:
        the parsed NNKP
    """
    comments = _core.extract_comments(string)
    parameters = _core.parse_parameters(_core.extract_parameters(string))
    blocks = _core.parse_blocks(_core.extract_blocks(string))

    parsed_nnkp = {
        'comments': comments,
        'parameters': parameters,
        'blocks': blocks,
        'direct_lattice': parse_direct_lattice(blocks['real_lattice']),
        'reciprocal_lattice': parse_reciprocal_lattice(blocks['recip_lattice']),
        'kpoints': parse_kpoints(blocks['kpoints']),
        'nnkpts': parse_nnkpts(blocks['nnkpts']),
        'exclude_bands': parse_exclude_bands(blocks['exclude_bands']),
    }
    if 'projections' in blocks:
        parsed_nnkp['projections'] = parse_projections(blocks['projections'])
    if 'spinor_projections' in blocks:
        parsed_nnkp['spinor_projections'] = parse_spinor_projections(blocks['spinor_projections'])

    return parsed_nnkp

wannier90io.parse_wout_iteration_info(stream: typing.TextIO) -> dict

Parse WOUT iteration info

Parameters:

Name Type Description Default
stream typing.TextIO

a file-like stream

required

Returns:

Type Description
dict

the parsed WOUT iteration info

Source code in wannier90io/_wout.py
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
def parse_wout_iteration_info(stream: typing.TextIO) -> dict:
    """
    Parse WOUT iteration info

    Arguments:
        stream: a file-like stream

    Returns:
        the parsed WOUT iteration info

    """
    convergence = []
    spread = []
    delta = []
    disentanglement = []
    for line in stream.readlines():
        match = patterns['convergence_line'].match(line)
        if match is not None:
            convergence += [{
                'iter': int(match.group('iter')),
                'delta': float(match.group('delta')),
                'gradient': float(match.group('gradient')),
                'spread': float(match.group('spread')),
                'time': float(match.group('time')),
            }]

        match = patterns['spread_line'].match(line)
        if match is not None:
            spread += [{
                'omega_d': float(match.group('omega_d')),
                'omega_od': float(match.group('omega_od')),
                'omega_tot': float(match.group('omega_tot')),
            }]

        match = patterns['delta_line'].match(line)
        if match is not None:
            delta += [{
                'omega_d': float(match.group('omega_d')),
                'omega_od': float(match.group('omega_od')),
                'omega_tot': float(match.group('omega_tot')),
            }]

        match = patterns['disentanglement_line'].match(line)
        if match is not None:
            disentanglement += [{
                'iter': int(match.group('iter')),
                'omega_i_i': float(match.group('omega_i_i')),
                'omega_i_f': float(match.group('omega_i_f')),
                'delta': float(match.group('delta')),
                'time': float(match.group('time')),
            }]

    return {
        'convergence': convergence,
        'spread': spread,
        'delta': delta,
        'disentanglement': disentanglement,
    }

wannier90io.read_amn(stream: typing.TextIO) -> np.ndarray

Read projections matrix

Parameters:

Name Type Description Default
stream typing.TextIO

a file-like stream

required

Returns:

Type Description
np.ndarray

projections matrix (Nk, Nb, Np)

Source code in wannier90io/_amn.py
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
def read_amn(stream: typing.TextIO) -> np.ndarray:
    """
    Read projections matrix

    Arguments:
        stream: a file-like stream

    Returns:
        projections matrix (Nk, Nb, Np)

    """
    stream.readline()

    [Nb, Nk, Np] = np.fromstring(stream.readline(), sep=' ', dtype=int)

    raw_data = np.loadtxt(stream).reshape((Nk, Np, Nb, 5))

    amn = np.transpose(raw_data[:, :, :, 3] + 1j*raw_data[:, :, :, 4], axes=(0, 2, 1))

    return amn

wannier90io.write_amn(stream: typing.TextIO, amn: np.ndarray, header: typing.Optional[str] = 'HEADER')

Write projections matrix

Parameters:

Name Type Description Default
stream typing.TextIO

a file-like stream

required
amn np.ndarray

projections matrix (Nk, Nb, Np)

required
header typing.Optional[str]

header

'HEADER'
Source code in wannier90io/_amn.py
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
def write_amn(stream: typing.TextIO, amn: np.ndarray, header: typing.Optional[str] = 'HEADER'):
    r"""
    Write projections matrix

    Arguments:
        stream: a file-like stream
        amn: projections matrix (Nk, Nb, Np)
        header: header

    """
    (Nk, Nb, Np) = amn.shape
    indices = np.mgrid[:Nb, :Np, :Nk].reshape((3, -1), order='F') + 1

    amn = np.transpose(amn, axes=(1, 2, 0)).flatten(order='F').view(float).reshape((-1, 2))
    data = np.column_stack((indices.T, amn))

    print(header, file=stream)
    print(f'{Nb:12d}{Nk:12d}{Np:12d}', file=stream)
    np.savetxt(stream, data, fmt='%5d%5d%5d%18.12f%18.12f')

wannier90io.read_chk(stream: typing.TextIO) -> dict

Read checkpoint

Parameters:

Name Type Description Default
stream typing.TextIO

a file-like stream

required

Returns:

Type Description
dict

dict

Source code in wannier90io/_chk.py
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
def read_chk(stream: typing.TextIO) -> dict:
    """
    Read checkpoint

    Arguments:
        stream: a file-like stream

    Returns:
        dict

    """
    chk = {}

    chk['header'] = stream.readline()
    chk['num_bands'] = Nb = int(stream.readline())
    chk['num_exclude_bands'] = int(stream.readline())
    if chk['num_exclude_bands'] > 0:
        chk['num_exclude_bands'] = np.fromstring(stream.readline(), dtype=int)
    chk['real_lattice'] = np.fromstring(stream.readline(), sep=' ', dtype=float).reshape((3, 3), order='F')
    chk['recip_lattice'] = np.fromstring(stream.readline(), sep=' ', dtype=float).reshape((3, 3), order='F')
    chk['num_kpts'] = Nk = int(stream.readline())
    chk['mp_grid'] = np.fromstring(stream.readline(), sep=' ', dtype=int)
    chk['kpt_latt'] = np.zeros((chk['num_kpts'], 3))
    for idx in range(chk['num_kpts']):
        chk['kpt_latt'][idx] = np.fromstring(stream.readline(), sep=' ', dtype=float)
    chk['nntot'] = Nn = int(stream.readline())
    chk['num_wann'] = Nw = int(stream.readline())
    chk['checkpoint'] = stream.readline()
    chk['have_disentangled'] = bool(int(stream.readline()))
    if chk['have_disentangled']:
        chk['omega_invariant'] = float(stream.readline())
        chk['lwindow'] = np.loadtxt(stream, max_rows=(Nk*Nb), dtype=bool).reshape((Nk, Nb))
        chk['nwindim'] = np.loadtxt(stream, max_rows=Nk, dtype=int)
        chk['u_matrix_opt'] = np.loadtxt(stream, max_rows=(Nk*Nw*Nb), dtype=float).view(complex).reshape((Nk, Nw, Nb))
    chk['u_matrix'] = np.loadtxt(stream, max_rows=(Nk*Nw*Nw), dtype=float).view(complex).reshape((Nw, Nw, Nk), order='F').transpose((2, 0, 1))
    chk['m_matrix'] = np.loadtxt(stream, max_rows=(Nk*Nn*Nw*Nw), dtype=float).view(complex).reshape((Nw, Nw, Nn, Nk), order='F').transpose((3, 2, 0, 1))

    return chk

wannier90io.read_u(stream: typing.TextIO) -> tuple[np.ndarray, np.ndarray]

Read unitary matrix file (seedname_u.mat) or the rectangular U_dis matrix file (seedname_u_dis.mat).

Note

for the _u.mat file, num_bands == num_wann.

Parameters:

Name Type Description Default
stream typing.TextIO

a file-like stream

required

Returns:

Type Description
np.ndarray

kpoint coordinates in fractional coordinates (num_kpts, 3)

np.ndarray

U matrix U(k) or U_dis(k) (num_kpts, num_bands, num_wann)

Source code in wannier90io/_u.py
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
def read_u(stream: typing.TextIO) -> tuple[np.ndarray, np.ndarray]:
    """
    Read unitary matrix file (seedname_u.mat) or the rectangular U_dis matrix
    file (seedname_u_dis.mat).

    Note:
        for the _u.mat file, num_bands == num_wann.

    Arguments:
        stream: a file-like stream

    Returns:
        kpoint coordinates in fractional coordinates (num_kpts, 3)
        U matrix U(k) or U_dis(k) (num_kpts, num_bands, num_wann)

    """
    stream.readline()   # header

    [nkpt, num_wann, num_bands] = np.fromstring(stream.readline(), sep=' ', dtype=int)
    u_matrices = np.zeros((nkpt, num_bands, num_wann), dtype=complex)
    kpoints = []

    for ikpt in range(nkpt):
        empty = stream.readline()   # header
        assert not empty.strip(), f"Expected empty line but found instead: '{empty}'"

        kpoint = np.fromstring(stream.readline(), sep=' ', dtype=float)
        assert len(kpoint) == 3
        kpoints.append(kpoint)
        u_matrices[ikpt, :, :] = np.loadtxt(stream, max_rows=(num_wann * num_bands)).view(complex).reshape((num_bands, num_wann), order='F')

    return np.array(kpoints), u_matrices

wannier90io.read_unk_formatted(stream: typing.TextIO) -> tuple[int, np.ndarray]

Read wavefunction files (UNKnnnnn.n files) in formatted format.

Note that the UNK files must have been created using the wvfn_formatted option set to True in the interface code (e.g. pw2wannier90.x for the Quantum ESPRESSO interface). Note that this is not the default, however for reading into an external code, this is recommended for portability.

Note

for now only works in the non-spinor case. Spinor case still to be implemented.

Parameters:

Name Type Description Default
stream typing.TextIO

a file-like stream

required

Returns:

Type Description
int

k-point index ik (integer)

np.ndarray

complex wavefunction (ngx, ngy, ngz, Nb)

Source code in wannier90io/_unk.py
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
def read_unk_formatted(stream: typing.TextIO) -> tuple[int, np.ndarray]:
    """
    Read wavefunction files (UNKnnnnn.n files) in formatted format.

    Note that the UNK files must have been created using the `wvfn_formatted`
    option set to True in the interface code (e.g. pw2wannier90.x for the
    Quantum ESPRESSO interface). Note that this is *not* the default, however
    for reading into an external code, this is recommended for portability.

    Note:
       for now only works in the non-spinor case.
       Spinor case still to be implemented.

    Arguments:
        stream: a file-like stream

    Returns:
        k-point index ik (integer)
        complex wavefunction (ngx, ngy, ngz, Nb)

    """
    [ngx, ngy, ngz, ik, nbnd] = np.fromstring(stream.readline(), sep=' ', dtype=int)

    wvfn = np.zeros((ngx, ngy, ngz, nbnd), dtype=complex)

    for ibnd in range(nbnd):
        wvfn[:, :, :, ibnd] = np.loadtxt(stream, max_rows=(ngx * ngy * ngz)).view(complex).reshape((ngx, ngy, ngz), order='F')

    return (ik, wvfn)

wannier90io.read_eig(stream: typing.TextIO) -> np.ndarray

Read eigenvalues matrix

Parameters:

Name Type Description Default
stream typing.TextIO

a file-like stream

required

Returns:

Type Description
np.ndarray

eigenvalues matrix (Nk, Nb)

Source code in wannier90io/_eig.py
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
def read_eig(stream: typing.TextIO) -> np.ndarray:
    """
    Read eigenvalues matrix

    Arguments:
        stream: a file-like stream

    Returns:
        eigenvalues matrix (Nk, Nb)

    """
    raw_data = np.loadtxt(stream)

    Nb = int(raw_data[-1, 0])
    Nk = int(raw_data[-1, 1])

    eig = raw_data[:, 2].reshape((Nk, Nb))

    return eig

wannier90io.write_eig(stream: typing.TextIO, eig: np.ndarray)

Write eigenvalues matrix

Parameters:

Name Type Description Default
stream typing.TextIO

a file-like stream

required
eig np.ndarray

eigenvalues matrix (Nk, Nb)

required
Source code in wannier90io/_eig.py
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
def write_eig(stream: typing.TextIO, eig: np.ndarray):
    r"""
    Write eigenvalues matrix

    Arguments:
        stream: a file-like stream
        eig: eigenvalues matrix (Nk, Nb)

    """
    (Nk, Nb) = eig.shape
    indices = np.mgrid[:Nb, :Nk].reshape((2, Nk*Nb), order='F') + 1

    data = np.column_stack((indices.T, eig.flatten()))

    np.savetxt(stream, data, fmt='%5d%5d%18.12f')

wannier90io.read_mmn(stream: typing.TextIO) -> tuple[np.ndarray, np.ndarray]

Read overlaps matrix

Parameters:

Name Type Description Default
stream typing.TextIO

a file-like stream

required

Returns:

Type Description
np.ndarray

overlaps matrix (Nk, Nn, Nb, Nb)

np.ndarray

nnkps (Nk, Nn, 5)

Source code in wannier90io/_mmn.py
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
def read_mmn(stream: typing.TextIO) -> tuple[np.ndarray, np.ndarray]:
    """
    Read overlaps matrix

    Arguments:
        stream: a file-like stream

    Returns:
        overlaps matrix (Nk, Nn, Nb, Nb)
        nnkps (Nk, Nn, 5)

    """
    stream.readline()   # header

    [Nb, Nk, Nn] = np.fromstring(stream.readline(), sep=' ', dtype=int)

    mmn = np.zeros((Nk, Nn, Nb, Nb), dtype=complex)
    nnkpts = np.zeros((Nk, Nn, 5), dtype=int)

    for (ik, ikb) in itertools.product(range(Nk), range(Nn)):
        nnkpts[ik, ikb] = np.fromstring(stream.readline(), sep=' ', dtype=int)
        mmn[ik, ikb] = np.loadtxt(stream, max_rows=(Nb*Nb)).view(complex).reshape((Nb, Nb), order='F')

    nnkpts[:, :, 0] -= 1
    nnkpts[:, :, 1] -= 1

    return (mmn, nnkpts)