Skip to content

Commit c4c980d

Browse files
Wentzellthe-hampel
authored andcommitted
Replace f2py Fortran module with pure Python for Elk PMAT.OUT reader
The elktools converter previously used an f2py-compiled Fortran module to read binary PMAT.OUT files. This has been replaced with a pure Python implementation that provides identical functionality without requiring a Fortran compiler at build time. Changes: - Add python/triqs_dft_tools/converters/elktools/getpmatelk.py Pure Python implementation to read Elk's Fortran unformatted binary files - Update python/triqs_dft_tools/converters/elk.py Modified import to use new Python module - Update python/triqs_dft_tools/converters/elktools/__init__.py Export getpmatelk function for direct import - Remove python/triqs_dft_tools/converters/elktools/elkwrappers/ Deleted f2py build infrastructure and Fortran source - Update CMakeLists.txt Removed elkwrappers subdirectory from build Benefits: - Eliminates Fortran compiler dependency for this component - Simplifies build process (no f2py compilation needed) - Easier to maintain and debug - Better error messages and handling - Resolves historical f2py compatibility issues on Intel/macOS Implementation details: - Correctly handles Fortran direct-access unformatted binary files - Properly reads column-major array ordering - Maintains exact compatibility with original Fortran version - All validation and error checking preserved Testing: - All 5 Elk converter tests pass - Output matches reference files exactly (verified with h5diff) - Tested with actual Elk PMAT.OUT files (41 states, 10 k-points) Cc: @AlynJ
1 parent 96fd210 commit c4c980d

File tree

7 files changed

+122
-110
lines changed

7 files changed

+122
-110
lines changed

CMakeLists.txt

Lines changed: 0 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -145,8 +145,6 @@ endif()
145145
# Python
146146
if(PythonSupport)
147147
add_subdirectory(python/${PROJECT_NAME})
148-
# elk binary i/o code for wrappers
149-
add_subdirectory(python/${PROJECT_NAME}/converters/elktools/elkwrappers)
150148
endif()
151149

152150
# Docs

python/triqs_dft_tools/converters/elk.py

Lines changed: 4 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -787,8 +787,8 @@ def convert_transport_input(self):
787787
# velocities_k: velocity (momentum) matrix elements between all bands in band_window_optics
788788
# and each k-point.
789789

790-
#load fortran wrapper module
791-
import triqs_dft_tools.converters.elktools.elkwrappers.getpmatelk as et
790+
#load pure Python module to read PMAT.OUT
791+
from triqs_dft_tools.converters.elktools import getpmatelk
792792
#elk velocities for all bands
793793
pmat=numpy.zeros([self.nstsv,self.nstsv,3],dtype=complex)
794794

@@ -802,10 +802,8 @@ def convert_transport_input(self):
802802
mpi.report("Reading PMAT.OUT")
803803
#read velocities for each k-point
804804
for ik in range(self.n_k):
805-
#need to use a fortran array for wrapper
806-
f_vkl = numpy.asfortranarray(self.vkl[ik,:])
807-
#read the ik velocity using the wrapper
808-
pmat[:,:,:]=et.getpmatelk(ik+1,self.nstsv,f_vkl)
805+
#read the ik velocity using pure Python
806+
pmat[:,:,:]=getpmatelk(ik+1, self.nstsv, self.vkl[ik,:])
809807
#loop over spin
810808
for isp in range(n_spin_blocks):
811809
#no. correlated bands at ik

python/triqs_dft_tools/converters/elktools/__init__.py

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -25,6 +25,7 @@
2525

2626
from .readElkfiles import readElkfiles
2727
from .elk_converter_tools import ElkConverterTools
28+
from .getpmatelk import getpmatelk
2829

29-
__all__ =['readElkfiles','ElkConverterTools']
30+
__all__ =['readElkfiles','ElkConverterTools','getpmatelk']
3031

python/triqs_dft_tools/converters/elktools/elkwrappers/CMakeLists.txt

Lines changed: 0 additions & 36 deletions
This file was deleted.

python/triqs_dft_tools/converters/elktools/elkwrappers/getpmatelk.f90

Lines changed: 0 additions & 65 deletions
This file was deleted.
Lines changed: 116 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,116 @@
1+
2+
##########################################################################
3+
#
4+
# TRIQS: a Toolbox for Research in Interacting Quantum Systems
5+
#
6+
# Copyright (C) 2025 by N. Wentzell
7+
#
8+
# TRIQS is free software: you can redistribute it and/or modify it under the
9+
# terms of the GNU General Public License as published by the Free Software
10+
# Foundation, either version 3 of the License, or (at your option) any later
11+
# version.
12+
#
13+
# TRIQS is distributed in the hope that it will be useful, but WITHOUT ANY
14+
# WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
15+
# FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
16+
# details.
17+
#
18+
# You should have received a copy of the GNU General Public License along with
19+
# TRIQS. If not, see <http://www.gnu.org/licenses/>.
20+
#
21+
##########################################################################
22+
23+
import numpy as np
24+
import os
25+
26+
def getpmatelk(ik, nstsv, vkl, filename='PMAT.OUT'):
27+
"""
28+
Pure Python implementation to read momentum matrix elements from Elk's PMAT.OUT file.
29+
30+
This replaces the previous f2py compiled Fortran module.
31+
32+
Parameters
33+
----------
34+
ik : int
35+
K-point index (1-based, as in Fortran)
36+
nstsv : int
37+
Number of states (second variation)
38+
vkl : array_like
39+
K-point lattice coordinates (3 floats)
40+
filename : str, optional
41+
Path to PMAT.OUT file (default: 'PMAT.OUT')
42+
43+
Returns
44+
-------
45+
pmat : ndarray
46+
Momentum matrix elements with shape (nstsv, nstsv, 3), dtype=complex128
47+
"""
48+
49+
if not os.path.exists(filename):
50+
raise IOError(f"Error(getpmat): unable to read from {filename}")
51+
52+
# Tolerance for k-point comparison (from Elk's epslat)
53+
epslat = 1.0e-6
54+
55+
# Data types matching Fortran:
56+
# real(8) -> float64 (8 bytes)
57+
# integer -> int32 (4 bytes, standard Fortran default integer)
58+
# Note: If Elk is compiled with -fdefault-integer-8, use int64
59+
# complex(8) -> complex128 (16 bytes: 8-byte real + 8-byte imag)
60+
61+
# Calculate record length in bytes
62+
# Record contains: vkl(3) + nstsv + pmat(nstsv, nstsv, 3)
63+
vkl_bytes = 3 * 8 # 3 real(8)
64+
nstsv_bytes = 4 # 1 integer
65+
pmat_bytes = nstsv * nstsv * 3 * 16 # nstsv*nstsv*3 complex(8)
66+
recl = vkl_bytes + nstsv_bytes + pmat_bytes
67+
68+
# Read the record at position ik (1-based indexing)
69+
# In direct access, record ik starts at byte position (ik-1) * recl
70+
try:
71+
with open(filename, 'rb') as f:
72+
# Seek to the correct record
73+
f.seek((ik - 1) * recl)
74+
75+
# Read vkl (3 float64)
76+
vkl_read = np.fromfile(f, dtype=np.float64, count=3)
77+
if vkl_read.size != 3:
78+
raise IOError(f"Error(getpmat): unable to read from {filename}")
79+
80+
# Read nstsv (1 int32)
81+
nstsv_array = np.fromfile(f, dtype=np.int32, count=1)
82+
if nstsv_array.size != 1:
83+
raise IOError(f"Error(getpmat): unable to read from {filename}")
84+
nstsv_read = nstsv_array[0]
85+
86+
# Read pmat (nstsv*nstsv*3 complex128)
87+
# Fortran stores arrays in column-major order
88+
# pmat(nstsv, nstsv, 3) in Fortran means first index varies fastest
89+
pmat_flat = np.fromfile(f, dtype=np.complex128, count=nstsv*nstsv*3)
90+
if pmat_flat.size != nstsv*nstsv*3:
91+
raise IOError(f"Error(getpmat): unable to read from {filename}")
92+
93+
# Reshape directly to (nstsv, nstsv, 3) with Fortran order
94+
pmat = pmat_flat.reshape((nstsv, nstsv, 3), order='F')
95+
except (OSError, IOError) as e:
96+
raise IOError(f"Error(getpmat): unable to read from {filename}") from e
97+
98+
# Validate k-point
99+
vkl = np.asarray(vkl)
100+
t1 = np.sum(np.abs(vkl - vkl_read))
101+
if t1 > epslat:
102+
raise ValueError(
103+
f"Error(getpmat): differing vectors for k-point {ik}\n"
104+
f" current : {vkl}\n"
105+
f" PMAT.OUT : {vkl_read}"
106+
)
107+
108+
# Validate nstsv
109+
if nstsv != nstsv_read:
110+
raise ValueError(
111+
f"Error(getpmat): differing nstsv for k-point {ik}\n"
112+
f" current : {nstsv}\n"
113+
f" PMAT.OUT : {nstsv_read}"
114+
)
115+
116+
return pmat
Binary file not shown.

0 commit comments

Comments
 (0)