-
Notifications
You must be signed in to change notification settings - Fork 59
/
pod.py
executable file
·168 lines (133 loc) · 4.96 KB
/
pod.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
"""
Module for Proper Orthogonal Decomposition (POD).
Three different methods can be employed: Truncated Singular Value Decomposition,
Truncated Randomized Singular Value Decomposition, Truncated Singular Value
Decomposition via correlation matrix.
"""
import numpy as np
from .reduction import Reduction
class POD(Reduction):
def __init__(self, method='svd', **kwargs):
"""
"""
available_methods = {
'svd': (self._svd, {
'rank': -1
}),
'randomized_svd': (self._rsvd, {
'rank': -1
}),
'correlation_matrix': (self._corrm, {
'rank': -1,
'save_memory': False
}),
}
self._modes = None
self._singular_values = None
method = available_methods.get(method)
if method is None:
raise RuntimeError(
"Invalid method for POD. Please chose one among {}".format(
', '.join(available_methods)))
self.__method, args = method
args.update(kwargs)
for hyperparam, value in args.items():
setattr(self, hyperparam, value)
@property
def modes(self):
"""
The POD modes.
:type: numpy.ndarray
"""
return self._modes
@property
def singular_values(self):
"""
The singular values
:type: numpy.ndarray
"""
return self._singular_values
def reduce(self, X):
"""
Reduces the parameter Space by using the specified reduction method (default svd).
:type: numpy.ndarray
"""
self._modes, self._singular_values = self.__method(X)
return self.modes.T.conj().dot(X)
def expand(self, X):
"""
Projects a reduced to full order solution.
:type: numpy.ndarray
"""
return self.modes.dot(X).T
def _truncation(self, X, s):
"""
Return the number of modes to select according to the `rank` parameter.
See POD.__init__ for further info.
:param numpy.ndarray X: the matrix to decompose.
:param numpy.ndarray s: the singular values of X.
:return: the number of modes
:rtype: int
"""
if self.rank == 0:
omega = lambda x: 0.56 * x**3 - 0.95 * x**2 + 1.82 * x + 1.43
beta = np.divide(*sorted(X.shape))
tau = np.median(s) * omega(beta)
rank = np.sum(s > tau)
elif self.rank > 0 and self.rank < 1:
cumulative_energy = np.cumsum(s**2 / (s**2).sum())
rank = np.searchsorted(cumulative_energy, self.rank) + 1
elif self.rank >= 1 and isinstance(self.rank, int):
rank = self.rank #min(self.rank, U.shape[1])
else:
rank = X.shape[1]
return rank
def _svd(self, X):
"""
Truncated Singular Value Decomposition.
:param numpy.ndarray X: the matrix to decompose.
:return: the truncated left-singular vectors matrix, the truncated
singular values array, the truncated right-singular vectors matrix.
:rtype: numpy.ndarray, numpy.ndarray, numpy.ndarray
References:
Gavish, Matan, and David L. Donoho, The optimal hard threshold for
singular values is, IEEE Transactions on Information Theory 60.8
(2014): 5040-5053.
"""
U, s = np.linalg.svd(X, full_matrices=False)[:2]
rank = self._truncation(X, s)
return U[:, :rank], s[:rank]
def _rsvd(self, X):
"""
Truncated randomized Singular Value Decomposition.
:param numpy.ndarray X: the matrix to decompose.
:return: the truncated left-singular vectors matrix, the truncated
singular values array, the truncated right-singular vectors matrix.
:rtype: numpy.ndarray, numpy.ndarray, numpy.ndarray
"""
P = np.random.rand(X.shape[1], X.shape[0])
Q = np.linalg.qr(X.dot(P))[0]
Y = Q.T.conj().dot(X)
Uy, s = np.linalg.svd(Y, full_matrices=False)[:2]
U = Q.dot(Uy)
rank = self._truncation(X, s)
return U[:, :rank], s[:rank]
def _corrm(self, X):
"""
Truncated Singular Value Decomposition. calculated with correlation matrix.
:param numpy.ndarray X: the matrix to decompose.
:return: the truncated left-singular vectors matrix, the truncated
singular values array, the truncated right-singular vectors matrix.
:rtype: numpy.ndarray, numpy.ndarray, numpy.ndarray
"""
if self.save_memory:
corr = np.empty(shape=(X.shape[1], X.shape[1]))
for i, i_snap in enumerate(X.T):
for j, k_snap in enumerate(X.T):
corr[i, j] = np.inner(i_snap, k_snap)
else:
corr = X.T.dot(X)
s, U = np.linalg.eig(corr)
U = X.dot(U) / np.sqrt(s)
rank = self._truncation(X, s)
return U[:, :rank], s[:rank]