/
random_matrix.py
202 lines (162 loc) · 6.4 KB
/
random_matrix.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
r"""Provide some random matrix ensembles for numpy.
The implemented ensembles are:
=========== ======================== ======================= ================== ===========
ensemble matrix class drawn from measure invariant under beta
=========== ======================== ======================= ================== ===========
GOE real, symmetric ``~ exp(-n/4 tr(H^2))`` orthogonal O 1
----------- ------------------------ ----------------------- ------------------ -----------
GUE hermitian ``~ exp(-n/2 tr(H^2))`` unitary U 2
----------- ------------------------ ----------------------- ------------------ -----------
CRE O(n) Haar orthogonal O /
----------- ------------------------ ----------------------- ------------------ -----------
COE U in U(n) with U = U^T Haar orthogonal O 1
----------- ------------------------ ----------------------- ------------------ -----------
CUE U(n) Haar unitary U 2
----------- ------------------------ ----------------------- ------------------ -----------
O_close_1 O(n) ? / /
----------- ------------------------ ----------------------- ------------------ -----------
U_close_1 U(n) ? / /
=========== ======================== ======================= ================== ===========
All functions in this module take a tuple ``(n, n)`` as first argument, such that
we can use the function :meth:`~tenpy.linalg.np_conserved.Array.from_func`
to generate a block diagonal :class:`~tenpy.linalg.np_conserved.Array` with the block from the
corresponding ensemble, for example::
npc.Array.from_func_square(GOE, [leg, leg.conj()])
"""
# Copyright 2018-2021 TeNPy Developers, GNU GPLv3
import numpy as np
__all__ = [
'box', 'standard_normal_complex', 'GOE', 'GUE', 'CRE', 'COE', 'CUE', 'O_close_1', 'U_close_1'
]
def box(size, W=1.):
"""return random number uniform in (-W, W]."""
return (0.5 - np.random.random(size)) * (2. * W)
def standard_normal_complex(size):
"""return ``(R + 1.j*I)`` for independent `R` and `I` from np.random.standard_normal."""
return np.random.standard_normal(size) + 1.j * np.random.standard_normal(size)
def GOE(size):
r"""Gaussian orthogonal ensemble (GOE).
Parameters
----------
size : tuple
``(n, n)``, where `n` is the dimension of the output matrix.
Returns
-------
H : ndarray
Real, symmetric numpy matrix drawn from the GOE, i.e.
:math:`p(H) = 1/Z exp(-n/4 tr(H^2))`
"""
A = np.random.standard_normal(size)
return (A + A.T) * 0.5
def GUE(size):
r"""Gaussian unitary ensemble (GUE).
Parameters
----------
size : tuple
``(n, n)``, where `n` is the dimension of the output matrix.
Returns
-------
H : ndarray
Hermitian (complex) numpy matrix drawn from the GUE, i.e.
:math:`p(H) = 1/Z exp(-n/4 tr(H^2))`.
"""
A = standard_normal_complex(size)
return (A + A.T.conj()) * 0.5
def CRE(size):
r"""Circular real ensemble (CRE).
Parameters
----------
size : tuple
``(n, n)``, where `n` is the dimension of the output matrix.
Returns
-------
U : ndarray
Orthogonal matrix drawn from the CRE (=Haar measure on O(n)).
"""
# almost same code as for CUE
n, m = size
assert n == m # ensure that `mode` in qr doesn't matter.
A = np.random.standard_normal(size)
Q, R = np.linalg.qr(A)
# Q-R is not unique; to make it unique ensure that the diagonal of R is positive
# Q' = Q*L; R' = L^{-1} *R, where L = diag(phase(diagonal(R)))
L = np.diagonal(R)
Q *= np.sign(L)
return Q
def COE(size):
r"""Circular orthogonal ensemble (COE).
Parameters
----------
size : tuple
``(n, n)``, where `n` is the dimension of the output matrix.
Returns
-------
U : ndarray
Unitary, symmetric (complex) matrix drawn from the COE (=Haar measure on this space).
"""
U = CUE(size)
return np.dot(U.T, U)
def CUE(size):
r"""Circular unitary ensemble (CUE).
Parameters
----------
size : tuple
``(n, n)``, where `n` is the dimension of the output matrix.
Returns
-------
U : ndarray
Unitary matrix drawn from the CUE (=Haar measure on U(n)).
"""
# almost same code as for CRE
n, m = size
assert n == m # ensure that `mode` in qr doesn't matter.
A = standard_normal_complex(size)
Q, R = np.linalg.qr(A)
# Q-R is not unique; to make it unique ensure that the diagonal of R is positive
# Q' = Q*L; R' = L^{-1} *R, where L = diag(phase(diagonal(R)))
L = np.diagonal(R).copy()
L[np.abs(L) < 1.e-15] = 1.
Q *= L / np.abs(L)
return Q
def O_close_1(size, a=0.01):
r"""return an random orthogonal matrix 'close' to the Identity.
Parameters
----------
size : tuple
``(n, n)``, where `n` is the dimension of the output matrix.
a : float
Parameter determining how close the result is on `O`;
:math:`\lim_{a \rightarrow 0} <|O-E|>_a = 0`` (where `E` is the identity).
Returns
-------
O : ndarray
Orthogonal matrix close to the identiy (for small `a`).
"""
n, m = size
assert n == m
A = GOE(size) / (2. * n)**0.5 # scale such that eigenvalues are in [-1, 1]
E = np.eye(size[0])
Q, R = np.linalg.qr(E + a * A)
L = np.diagonal(R) # make QR decomposition unique & ensure Q is close to one for small `a`
Q *= np.sign(L)
return Q
def U_close_1(size, a=0.01):
r"""return an random orthogonal matrix 'close' to the identity.
Parameters
----------
size : tuple
``(n, n)``, where `n` is the dimension of the output matrix.
a : float
Parameter determining how close the result is to the identity.
:math:`\lim_{a \rightarrow 0} <|O-E|>_a = 0`` (where `E` is the identity).
Returns
-------
U : ndarray
Unitary matrix close to the identiy (for small `a`).
Eigenvalues are chosen i.i.d. as ``exp(1.j*a*x)`` with `x` uniform in [-1, 1].
"""
n, m = size
assert n == m
U = CUE(size) # random unitary
E = np.exp(1.j * a * (np.random.rand(n) * 2. - 1.))
return np.dot(U * E, U.T.conj())