-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathmk_preset.py
340 lines (306 loc) · 11 KB
/
mk_preset.py
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
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
import sparse_ir
import xprec # Make sure xprec is install for good accuracy
import numpy as np
# max line continuation is 255 for Fortran >= 2003,
# 39 for Fortran95.
# For safety, we stick to Fortran95.
MAX_CONTINUATION_LINE = 39
NLINE = MAX_CONTINUATION_LINE - 2
class BasisInfo:
def __init__(self, nlambda: int, ndigit: int) -> None:
beta = 1.0
lambda_ = 10.0 ** nlambda
wmax = lambda_/beta
eps = 1/10.0 ** ndigit
#kernel = sparse_ir.LogisticKernel(lambda_)
#sve_result = sparse_ir.compute_sve(kernel, eps)
#self.basis_f = sparse_ir.DimensionlessBasis("F", lambda_, eps, kernel=kernel, sve_result=sve_result)
#self.basis_b = sparse_ir.DimensionlessBasis("B", lambda_, eps, kernel=kernel, sve_result=sve_result)
basis = sparse_ir.FiniteTempBasisSet(1.0, lambda_, eps)
self.basis_f = basis.basis_f
self.basis_b = basis.basis_b
self.s = self.basis_f.s / np.sqrt(0.5 * lambda_)
self.size = self.s.size
# Sampling points:
# Add tau = 0, beta to sampling points (for Matsubara summation)
self.tau = np.unique(np.hstack([0, self.basis_f.default_tau_sampling_points(), beta]))
self.x = 2 * self.tau/beta - 1
self.freq_f = self.basis_f.default_matsubara_sampling_points()
self.freq_b = self.basis_b.default_matsubara_sampling_points()
self.omega = np.unique(np.hstack([-wmax, self.basis_f.default_omega_sampling_points(), wmax]))
self.y = self.omega/wmax
self.ntau = self.tau.size
self.nfreq_f = self.freq_f.size
self.nfreq_b = self.freq_b.size
self.nomega = self.omega.size
self.ntau_reduced = self.ntau //2 + 1
self.nfreq_f_reduced = self.nfreq_f //2 + 1
self.nfreq_b_reduced = self.nfreq_b //2 + 1
self.nomega_reduced = self.nomega //2 + 1
# Transformation matrix
self.u = np.sqrt(0.5 * beta) * self.basis_f.u(self.tau).T
self.uhat_f = (1/np.sqrt(beta)) * self.basis_f.uhat(self.freq_f).T
self.uhat_b = (1/np.sqrt(beta)) * self.basis_b.uhat(self.freq_b).T
self.v = np.sqrt(wmax) * self.basis_f.v(self.omega).T
def _str(value):
if isinstance(value, float):
return f"{value:.16e}".replace('e', 'd')
elif type(value) in [int, np.int64]:
return str(value)
else:
raise RuntimeError("Invalid type: " + str(type(value)))
def run(nlambda_list, ndigit_list):
bases = {}
for nlambda in nlambda_list:
for ndigit in ndigit_list:
bases[(nlambda, ndigit)] = BasisInfo(nlambda, ndigit)
nlambda_str = ' '.join(map(str, nlambda_list))
ndigit_str = ' '.join(map(str, ndigit_list))
print(
f"""\
!! Precomputed data of IR for a parameter matrix of
!! nlambda = {nlambda_str} and ndigit = {ndigit_str}.
!! This file was automatically generated:
!! python3 mk_preset.py --nlambda {nlambda_str} --ndigit {ndigit_str} > sparse_ir_preset.f90
!! Do not edit this file manually.
!!
!-----------------------------------------------------------------------
MODULE sparse_ir_preset
!-----------------------------------------------------------------------
USE sparse_ir
!
IMPLICIT NONE
!
PRIVATE
!
INTEGER, PARAMETER :: DP = KIND(0d0)
!
PUBLIC :: mk_ir_preset
!
""", end="")
for nlambda in nlambda_list:
for ndigit in ndigit_list:
b = bases[(nlambda, ndigit)]
sig = f"nlambda{nlambda}_ndigit{ndigit}"
print(
f"""\
REAL(KIND = DP) :: s_{sig}({b.size})
REAL(KIND = DP) :: tau_{sig}({b.ntau})
INTEGER :: freq_f_{sig}({b.nfreq_f})
INTEGER :: freq_b_{sig}({b.nfreq_b})
REAL(KIND = DP) :: omega_{sig}({b.nomega})
REAL(KIND = DP) :: u_r_{sig}({b.ntau_reduced} * {b.size})
REAL(KIND = DP) :: uhat_f_r_{sig}({b.nfreq_f_reduced} * {b.size})
REAL(KIND = DP) :: uhat_b_r_{sig}({b.nfreq_b_reduced} * {b.size})
REAL(KIND = DP) :: v_r_{sig}({b.nomega_reduced} * {b.size})
""", end="")
print(
"""\
!
CONTAINS
!
!-----------------------------------------------------------------------
FUNCTION mk_ir_preset(nlambda, ndigit, beta, positive_only) result(obj)
!-----------------------------------------------------------------------
!
INTEGER, INTENT(IN) :: nlambda, ndigit
REAL(KIND = DP), INTENT(IN) :: beta
LOGICAL, INTENT(IN), OPTIONAL :: positive_only
TYPE(IR) :: obj
!
""", end="")
for nlambda in nlambda_list:
for ndigit in ndigit_list:
print(
f"""\
IF (nlambda == {nlambda} .and. ndigit == {ndigit}) THEN
IF ((.NOT. PRESENT(positive_only))) THEN
obj = mk_nlambda{nlambda}_ndigit{ndigit}(beta)
ELSE
obj = mk_nlambda{nlambda}_ndigit{ndigit}(beta, positive_only)
ENDIF
RETURN
ENDIF
!
""", end="")
print(
"""\
STOP "Invalid parameters"
!
!-----------------------------------------------------------------------
END FUNCTION mk_ir_preset
!-----------------------------------------------------------------------
!
""", end="")
for nlambda in nlambda_list:
for ndigit in ndigit_list:
print_data(nlambda, ndigit, bases[(nlambda, ndigit)])
print(
"""\
!-----------------------------------------------------------------------
END MODULE sparse_ir_preset
!-----------------------------------------------------------------------
""", end="")
def print_data(nlambda, ndigit, b):
sig = f"nlambda{nlambda}_ndigit{ndigit}"
print(
f"""\
!-----------------------------------------------------------------------
FUNCTION mk_nlambda{nlambda}_ndigit{ndigit}(beta, positive_only) result(obj)
!-----------------------------------------------------------------------
!
REAL(KIND = DP), INTENT(IN) :: beta
LOGICAL, INTENT(IN), OPTIONAL :: positive_only
TYPE(IR) :: obj
REAL(KIND = DP), ALLOCATABLE :: u(:, :), v(:, :), dlr(:, :)
COMPLEX(KIND = DP), ALLOCATABLE :: uhat_f(:, :), uhat_b(:, :)
INTEGER, PARAMETER :: size = {b.size}, ntau = {b.ntau}, nfreq_f = {b.nfreq_f}, nfreq_b = {b.nfreq_b}, nomega = {b.nomega}
INTEGER, PARAMETER :: nlambda = {nlambda}, ndigit = {ndigit}
INTEGER, PARAMETER :: ntau_reduced = ntau/2+1
INTEGER, PARAMETER :: nfreq_f_reduced = nfreq_f/2+1, nfreq_b_reduced = nfreq_b/2+1
INTEGER, PARAMETER :: nomega_reduced = nomega/2+1
REAL(KIND = DP), PARAMETER :: lambda = 1.d{nlambda}, eps = 1.d-{ndigit}
!
INTEGER :: itau, l, ifreq, iomega
!
""", end="")
for varname in ["s", "tau", "freq_f", "freq_b", "u_r", "uhat_f_r", "uhat_b_r", "omega", "v_r"]:
print(2*" " + f"CALL init_{varname}_{sig}()")
print(
f"""\
ALLOCATE(u(ntau, size))
ALLOCATE(uhat_f(nfreq_f, size))
ALLOCATE(uhat_b(nfreq_b, size))
ALLOCATE(v(nomega, size))
ALLOCATE(dlr(nomega, size))
!
! Use the fact U_l(tau) is even/odd for even/odd l-1.
DO l = 1, size
DO itau = 1, ntau_reduced
u(itau, l) = u_r_{sig}(itau + {b.ntau_reduced}*(l-1))
u(ntau-itau+1, l) = (-1)**(l-1) * u_r_{sig}(itau + {b.ntau_reduced}*(l-1))
ENDDO
ENDDO
!
! Use the fact U^F_l(iv) is pure imaginary/real for even/odd l-1.
DO l = 1, size, 2
DO ifreq = 1, nfreq_f_reduced
uhat_f(ifreq, l) = CMPLX(0.0d0, uhat_f_r_{sig}(ifreq + {b.nfreq_f_reduced}*(l-1)), KIND = DP)
ENDDO
ENDDO
DO l = 2, size, 2
DO ifreq = 1, nfreq_f_reduced
uhat_f(ifreq, l) = CMPLX(uhat_f_r_{sig}(ifreq + {b.nfreq_f_reduced}*(l-1)), 0.0, KIND = DP)
ENDDO
ENDDO
DO l = 1, size
DO ifreq = 1, nfreq_f
uhat_f(nfreq_f-ifreq+1, l) = CONJG(uhat_f(ifreq, l))
ENDDO
ENDDO
!
! Use the fact U^B_l(iv) is pure real/imaginary for even/odd l-1
DO l = 1, size, 2
DO ifreq = 1, nfreq_b_reduced
uhat_b(ifreq, l) = CMPLX(uhat_b_r_{sig}(ifreq + {b.nfreq_b_reduced}*(l-1)), 0.0d0, KIND = DP)
ENDDO
ENDDO
DO l = 2, size, 2
DO ifreq = 1, nfreq_b_reduced
uhat_b(ifreq, l) = CMPLX(0.0d0, uhat_b_r_{sig}(ifreq + {b.nfreq_b_reduced}*(l-1)), KIND = DP)
ENDDO
ENDDO
DO l = 1, size
DO ifreq = 1, nfreq_b
uhat_b(nfreq_b-ifreq+1, l) = CONJG(uhat_b(ifreq, l))
ENDDO
ENDDO
!
! Use the fact V_l(omega) is even/odd for even/odd l-1.
DO l = 1, size
DO iomega = 1, nomega_reduced
v(iomega, l) = v_r_{sig}(iomega + {b.nomega_reduced}*(l-1))
v(nomega-iomega+1, l) = (-1)**(l-1) * v_r_{sig}(iomega + {b.nomega_reduced}*(l-1))
ENDDO
DO iomega = 1, nomega
dlr(iomega, l) = - s_{sig}(l) * v(iomega, l)
ENDDO
ENDDO
!
IF ((.NOT. PRESENT(positive_only))) THEN
CALL init_ir(obj, beta, lambda, eps,&
s_{sig}, tau_{sig},&
freq_f_{sig}, freq_b_{sig},&
u, uhat_f, uhat_b, omega_{sig},&
v, dlr, 1d-20)
ELSE
CALL init_ir(obj, beta, lambda, eps,&
s_{sig}, tau_{sig},&
freq_f_{sig}, freq_b_{sig},&
u, uhat_f, uhat_b, omega_{sig},&
v, dlr, 1d-20, positive_only)
ENDIF
!
DEALLOCATE(u, uhat_f, uhat_b, v, dlr)
!
!-----------------------------------------------------------------------
END FUNCTION mk_nlambda{nlambda}_ndigit{ndigit}
!-----------------------------------------------------------------------
!
""", end="")
print_vector_data(b.s, f"s_{sig}")
print_vector_data(b.x, f"tau_{sig}")
print_vector_data(b.freq_f, f"freq_f_{sig}")
print_vector_data(b.freq_b, f"freq_b_{sig}")
print_vector_data(b.y, f"omega_{sig}")
print_vector_data(b.u[0:b.ntau_reduced,:], f"u_r_{sig}")
print_vector_data((b.uhat_f.real + b.uhat_f.imag)[0:b.nfreq_f_reduced,:], f"uhat_f_r_{sig}")
print_vector_data((b.uhat_b.real + b.uhat_b.imag)[0:b.nfreq_b_reduced,:], f"uhat_b_r_{sig}")
print_vector_data(b.v[0:b.nomega_reduced,:], f"v_r_{sig}")
def _array_to_strings(arr, num_elem_str = 3):
""" Convert array of float or int to a list of strings """
start = 0
res = []
while start < arr.size:
end = min(start + num_elem_str, arr.size)
res.append(", ".join(map(_str, arr[start:end])))
start = end
return res
def print_vector_data(vec, var_name):
""" Print initializer of array data """
vec = np.asfortranarray(vec)
vec = vec.T.ravel()
print(
f"""\
!-----------------------------------------------------------------------
SUBROUTINE init_{var_name}()
!-----------------------------------------------------------------------
!
""", end="")
start = 0
while start < vec.size:
end = min(start + NLINE, vec.size)
sub_vec = vec[start:end]
print(2*" " + f"{var_name}({start+1}:{end}) = (/ &")
end_str = "&"
print(", &\n".join(_array_to_strings(sub_vec)) + end_str)
print(2*" " + "/)")
start = end
print(
f"""\
!
!-----------------------------------------------------------------------
END SUBROUTINE init_{var_name}
!-----------------------------------------------------------------------
!
""", end="")
if __name__ == '__main__':
import argparse
p = argparse.ArgumentParser()
# accept two lists of arguments
p.add_argument('--nlambda', nargs="+", type=int)
p.add_argument('--ndigit', nargs="+", type=int)
args = p.parse_args()
nlambda_list = np.array(args.nlambda, dtype=int)
ndigit_list = np.array(args.ndigit, dtype=int)
run(nlambda_list, ndigit_list)