forked from wannier-developers/wannier90
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathc_interface.F90
305 lines (276 loc) · 11.3 KB
/
c_interface.F90
1
2
3
4
5
6
7
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
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
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
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
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
module w90_library_c
!Fortran 2018: Assumed-length character dummy argument ‘keyword’ at (1) of procedure ‘cset_option_int’ with BIND(C) attribute
use iso_c_binding
use w90_library
implicit none
public
type, bind(c) :: w90_data
type(c_ptr) :: caddr
end type
contains
subroutine w90_create(w90_obj) bind(c)
!! return a c-pointer to a instance of the wannier90 library data structure
type(lib_common_type), pointer :: common_data
type(w90_data) :: w90_obj
if (c_associated(w90_obj%caddr)) return
allocate (common_data)
w90_obj%caddr = c_loc(common_data)
end subroutine
subroutine w90_delete(w90_obj) bind(c)
!! deallocates/clears a c-pointer to a instance of the wannier90 library data structure
implicit none
type(w90_data) :: w90_obj
type(lib_common_type), pointer :: w90_fptr
if (.not. c_associated(w90_obj%caddr)) return
call c_f_pointer(w90_obj%caddr, w90_fptr)
deallocate (w90_fptr)
w90_obj%caddr = C_NULL_PTR
end subroutine
subroutine w90_disentangle_c(w90_obj, ierr) bind(c)
implicit none
type(w90_data), value :: w90_obj
type(lib_common_type), pointer :: w90_fptr
integer(kind=c_int) :: istdout, istderr, ierr
call w90_get_fortran_stderr(istderr)
call w90_get_fortran_stdout(istdout)
call c_f_pointer(w90_obj%caddr, w90_fptr)
call w90_disentangle(w90_fptr, istdout, istderr, ierr)
end subroutine
subroutine w90_get_centres_c(w90_obj, centres) bind(c)
! returns the centres of calulated mlwfs
implicit none
type(w90_data), value :: w90_obj
type(c_ptr), value :: centres
type(lib_common_type), pointer :: w90_fptr
real(kind=8), pointer :: fcentres(:, :)
call c_f_pointer(w90_obj%caddr, w90_fptr)
call c_f_pointer(centres, fcentres, [3, w90_fptr%num_wann])
call w90_get_centres(w90_fptr, fcentres)
end subroutine
subroutine w90_get_gkpb_c(w90_obj, gkpb) bind(c)
! return the g-offset of adjacent k-points in finite difference scheme
implicit none
type(w90_data), value :: w90_obj
type(c_ptr), value :: gkpb
type(lib_common_type), pointer :: w90_fptr
integer(kind=c_int), pointer :: nfptr(:, :, :)
integer(kind=c_int) :: istderr, istdout, ierr
call w90_get_fortran_stderr(istderr)
call w90_get_fortran_stdout(istdout)
call c_f_pointer(w90_obj%caddr, w90_fptr)
call c_f_pointer(gkpb, nfptr, [3, w90_fptr%num_kpts, w90_fptr%kmesh_info%nntot])
call w90_get_gkpb(w90_fptr, nfptr, istdout, istderr, ierr)
end subroutine
subroutine w90_get_nn_c(w90_obj, n) bind(c)
! return the number of adjacent k-points in finite difference scheme
implicit none
type(w90_data), value :: w90_obj
type(c_ptr), value :: n
type(lib_common_type), pointer :: w90_fptr
integer(kind=c_int), pointer :: ndat
integer(kind=c_int) :: istderr, istdout, ierr
call w90_get_fortran_stderr(istderr)
call w90_get_fortran_stdout(istdout)
call c_f_pointer(w90_obj%caddr, w90_fptr)
call c_f_pointer(n, ndat)
call w90_get_nn(w90_fptr, ndat, istdout, istderr, ierr)
end subroutine
subroutine w90_get_nnkp_c(w90_obj, nnkp) bind(c)
! return the indexing of adjacent k-points in finite difference scheme
implicit none
type(w90_data), value :: w90_obj
type(c_ptr), value :: nnkp
type(lib_common_type), pointer :: w90_fptr
integer(kind=c_int), pointer :: nfptr(:, :)
integer(kind=c_int) :: istderr, istdout, ierr
call w90_get_fortran_stderr(istderr)
call w90_get_fortran_stdout(istdout)
call c_f_pointer(w90_obj%caddr, w90_fptr)
call c_f_pointer(nnkp, nfptr, [w90_fptr%num_kpts, w90_fptr%kmesh_info%nntot])
call w90_get_nnkp(w90_fptr, nfptr, istdout, istderr, ierr)
end subroutine
subroutine w90_get_spreads_c(w90_obj, spreads) bind(c)
! returns the spreads of calulated mlwfs
implicit none
type(w90_data), value :: w90_obj
type(c_ptr), value :: spreads
type(lib_common_type), pointer :: w90_fptr
real(kind=8), pointer :: fspreads(:)
call c_f_pointer(w90_obj%caddr, w90_fptr)
call c_f_pointer(spreads, fspreads, [w90_fptr%num_wann])
call w90_get_spreads(w90_fptr, fspreads)
end subroutine
subroutine w90_input_setopt_c(w90_obj, seedname, ierr) bind(c)
! specify parameters through the library interface
implicit none
type(w90_data), value :: w90_obj
character(*, kind=c_char) :: seedname
type(lib_common_type), pointer :: w90_fptr
integer(kind=c_int) :: istderr, istdout, ierr
call w90_get_fortran_stderr(istderr)
call w90_get_fortran_stdout(istdout)
call c_f_pointer(w90_obj%caddr, w90_fptr)
call w90_input_setopt(w90_fptr, seedname, istdout, istderr, ierr)
end subroutine
subroutine w90_input_reader_c(w90_obj, ierr) bind(c)
! read (optional) parameters from .win file
implicit none
type(w90_data), value :: w90_obj
type(lib_common_type), pointer :: w90_fptr
integer(kind=c_int) :: istderr, istdout, ierr
call w90_get_fortran_stderr(istderr)
call w90_get_fortran_stdout(istdout)
call c_f_pointer(w90_obj%caddr, w90_fptr)
call w90_input_reader(w90_fptr, istdout, istderr, ierr)
end subroutine
subroutine w90_project_overlap_c(w90_obj, ierr) bind(c)
implicit none
type(w90_data), value :: w90_obj
type(lib_common_type), pointer :: w90_fptr
integer(kind=c_int) :: istdout, istderr, ierr
call w90_get_fortran_stderr(istderr)
call w90_get_fortran_stdout(istdout)
call c_f_pointer(w90_obj%caddr, w90_fptr)
call w90_project_overlap(w90_fptr, istdout, istderr, ierr)
end subroutine
subroutine w90_set_eigval_c(w90_obj, eigval_cptr) bind(c)
! copy a pointer to eigenvalue data
implicit none
type(w90_data), value :: w90_obj
type(c_ptr), value :: eigval_cptr
real(8), pointer :: eigval_fptr(:, :)
type(lib_common_type), pointer :: w90_fptr
call c_f_pointer(w90_obj%caddr, w90_fptr)
call c_f_pointer(eigval_cptr, eigval_fptr, [w90_fptr%num_bands, w90_fptr%num_kpts])
call w90_set_eigval(w90_fptr, eigval_fptr)
end subroutine
subroutine w90_set_m_local_c(w90_obj, m_cptr) bind(c)
! copy a pointer to m-matrix data
implicit none
type(w90_data), value :: w90_obj
type(c_ptr), value :: m_cptr
complex(8), pointer :: fptr(:, :, :, :)
type(lib_common_type), pointer :: w90_fptr
call c_f_pointer(w90_obj%caddr, w90_fptr)
call c_f_pointer(m_cptr, fptr, [w90_fptr%num_bands, w90_fptr%num_bands, w90_fptr%kmesh_info%nntot, w90_fptr%num_kpts])
call w90_set_m_local(w90_fptr, fptr)
end subroutine
subroutine w90_set_u_matrix_c(w90_obj, a_cptr) bind(c)
! copy pointer to u-matrix
implicit none
type(w90_data), value :: w90_obj
type(c_ptr), value :: a_cptr
complex(8), pointer :: a_fptr(:, :, :)
type(lib_common_type), pointer :: w90_fptr
call c_f_pointer(w90_obj%caddr, w90_fptr)
call c_f_pointer(a_cptr, a_fptr, [w90_fptr%num_wann, w90_fptr%num_wann, w90_fptr%num_kpts]) ! these are reversed wrt c
call w90_set_u_matrix(w90_fptr, a_fptr)
end subroutine
subroutine w90_set_u_opt_c(w90_obj, a_cptr) bind(c)
! copy pointer to u-matrix (also used for initial projections)
implicit none
type(w90_data), value :: w90_obj
type(c_ptr), value :: a_cptr
complex(kind=8), pointer :: a_fptr(:, :, :)
type(lib_common_type), pointer :: w90_fptr
call c_f_pointer(w90_obj%caddr, w90_fptr)
call c_f_pointer(a_cptr, a_fptr, [w90_fptr%num_bands, w90_fptr%num_wann, w90_fptr%num_kpts]) ! these are reversed wrt c
call w90_set_u_opt(w90_fptr, a_fptr)
end subroutine
subroutine w90_wannierise_c(w90_obj, ierr) bind(c)
implicit none
type(w90_data), value :: w90_obj
type(lib_common_type), pointer :: w90_fptr
integer(kind=c_int) :: istdout, istderr, ierr
call w90_get_fortran_stderr(istderr)
call w90_get_fortran_stdout(istdout)
call c_f_pointer(w90_obj%caddr, w90_fptr)
call w90_wannierise(w90_fptr, istdout, istderr, ierr)
end subroutine
subroutine w90_set_option_double_c(w90_obj, keyword, cdble) bind(c)
implicit none
type(w90_data), value :: w90_obj
character(*, kind=c_char) :: keyword
real(kind=c_double), value :: cdble
type(lib_common_type), pointer :: w90_fptr
call c_f_pointer(w90_obj%caddr, w90_fptr)
call w90_set_option(w90_fptr, keyword, cdble)
end subroutine
subroutine w90_set_option_double1d_c(w90_obj, keyword, arg_cptr, x) bind(c)
implicit none
type(w90_data), value :: w90_obj
character(*, kind=c_char) :: keyword
type(c_ptr), value :: arg_cptr
integer(kind=c_int), value :: x
type(lib_common_type), pointer :: w90_fptr
real(kind=8), pointer :: fptr(:)
call c_f_pointer(arg_cptr, fptr, [x])
call c_f_pointer(w90_obj%caddr, w90_fptr)
call w90_set_option(w90_fptr, keyword, fptr)
end subroutine
subroutine w90_set_option_double2d_c(w90_obj, keyword, arg_cptr, x, y) bind(c)
implicit none
type(w90_data), value :: w90_obj
character(*, kind=c_char) :: keyword
type(c_ptr), value :: arg_cptr
integer(kind=c_int), value :: x, y
type(lib_common_type), pointer :: w90_fptr
real(kind=8), pointer :: fptr(:, :)
call c_f_pointer(arg_cptr, fptr, [y, x]) ! these are reversed wrt c
call c_f_pointer(w90_obj%caddr, w90_fptr)
call w90_set_option(w90_fptr, keyword, fptr)
end subroutine
subroutine w90_set_option_int_c(w90_obj, keyword, cint) bind(c)
implicit none
type(w90_data), value :: w90_obj
character(*, kind=c_char) :: keyword
integer(kind=c_int), value :: cint
type(lib_common_type), pointer :: w90_fptr
call c_f_pointer(w90_obj%caddr, w90_fptr)
call w90_set_option(w90_fptr, keyword, cint)
end subroutine
subroutine w90_set_option_int1d_c(w90_obj, keyword, arg_cptr, x) bind(c)
implicit none
type(w90_data), value :: w90_obj
type(c_ptr), value :: arg_cptr
character(*, kind=c_char) :: keyword
integer(kind=c_int), value :: x
integer, pointer :: fptr(:)
type(lib_common_type), pointer :: w90_fptr
call c_f_pointer(w90_obj%caddr, w90_fptr)
call c_f_pointer(arg_cptr, fptr, [x])
call w90_set_option(w90_fptr, keyword, fptr)
end subroutine
subroutine w90_set_option_int2d_c(w90_obj, keyword, arg_cptr, x, y) bind(c)
implicit none
type(w90_data), value :: w90_obj
type(c_ptr), value :: arg_cptr
character(*, kind=c_char) :: keyword
integer(kind=c_int), value :: x, y
integer(kind=c_int), pointer :: fptr(:, :)
type(lib_common_type), pointer :: w90_fptr
call c_f_pointer(w90_obj%caddr, w90_fptr)
call c_f_pointer(arg_cptr, fptr, [y, x]) ! these are reversed wrt c
call w90_set_option(w90_fptr, keyword, fptr)
end subroutine
subroutine w90_set_option_logical_c(w90_obj, keyword, bool) bind(c)
implicit none
type(w90_data), value :: w90_obj
character(*, kind=c_char) :: keyword
logical(kind=c_bool), value :: bool
logical :: fbool
type(lib_common_type), pointer :: w90_fptr
fbool = bool
call c_f_pointer(w90_obj%caddr, w90_fptr)
call w90_set_option(w90_fptr, keyword, fbool)
end subroutine
subroutine w90_set_option_text_c(w90_obj, keyword, text) bind(c)
implicit none
type(w90_data), value :: w90_obj
character(*, kind=c_char) :: keyword
character(*, kind=c_char) :: text
type(lib_common_type), pointer :: w90_fptr
call c_f_pointer(w90_obj%caddr, w90_fptr)
call w90_set_option(w90_fptr, keyword, text)
end subroutine
end module w90_library_c