-
Notifications
You must be signed in to change notification settings - Fork 1
/
varpro.py
382 lines (345 loc) · 15.3 KB
/
varpro.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
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
import numpy as np
from scipy.sparse import spdiags
from scipy.optimize import least_squares
from scipy import linalg
import logging
logger = logging.getLogger(__name__)
class RankError(Exception):
'''Raised when linear parameters are ill-determined.'''
def __init__(self, rank, num_c, linear_params, nonlin_params):
message = """
The linear parameters are not well-determined.
The rank of the matrix in the subproblem is %d,
which is less than the number of linear parameters (%d).
Linear parameters:
%s
Nonlinear parameters:
%s""" \
% (rank, num_c, " ".join(map(str, linear_params)),
" ".join(map(str, nonlin_params)))
super().__init__(message)
class MaxIterationError(Exception):
'''Raised when maximum number of iterations is reached
during nonlinear solve.'''
def __init__(self, linear_params, nonlin_params):
message = """
Max iterations reached without convergence.
Linear parameters:
%s
Nonlinear parameters:
%s""" \
% (" ".join(map(str, linear_params)),
" ".join(map(str, nonlin_params)))
super().__init__(message)
def varpro(t, y, w, alpha, n, ada,
bounds = None, max_nfev = None, **kwargs):
#Solve a separable nonlinear least squares problem.
# This is a slightly simplified Python translation of the Matlab code by
# Dianne P. O'Leary and Bert W. Rust. The original code is documented in
# Computational Optimization and Applications, 54, 579 (2013). The Matlab
# code itself with extensive comments and references is available in the
# supplementary material to the online version of this paper at
# doi:10.1007/s10589-012-9492-9
#
# Given a set of m observations y(0),...,y(m-1) at "times" t(0),...,t(m-1),
# this program computes a weighted least squares fit using the model
#
# eta(alpha,c,t) =
# c_0 * phi_0 (alpha,t) + ... + c_(n-1) * phi_(n-1) (alpha,t)
# (possibly with an extra term + phi_n (alpha,t) ).
#
# This program determines optimal values of the q nonlinear parameters
# alpha and the n linear parameters c, given observations y at m
# different values of the "time" t and given evaluation of phi and
# derivatives of phi.
#
# On Input:
#
# t 1-d array containing the m "times" (independent variables).
# y 1-d array containing the m observations.
# w 1-d array containing the m weights used in the least squares
# fit. We minimize the norm of the weighted residual
# vector r, where, for i=0:m-1,
#
# r(i) = w(i) * (y(i) - eta(alpha, c, t(i))).
#
# Therefore, w(i) should be set to 1 divided by
# the standard deviation in the measurement y(i).
# If this number is unknown, set w(i) = 1.
# alpha 1-d array of length q with initial estimates of the parameters alpha.
# n number of linear parameters c
# ada a function described below.
# bounds If supplied, must be a tuple of both lower and upper
# (Optional) bounds for all q of the parameters alpha. See documentation
# of scipy.optimize.least_squares for examples.
# **kwargs Any keyword argument taken by scipy.optimize.least_squares
# (Optional) other than bounds can be passed here. See
# scipy.optimize.least_squares documentation for possible
# arguments.
#
# On Output:
#
# alpha length q Estimates of the nonlinear parameters.
# c length n Estimates of the linear parameters.
# wresid length m Weighted residual vector, with i-th component
# w(i) * (y(i) - eta(alpha, c, t(i))).
# wresid_norm Norm of wresid.
# y_est length m The model estimates = eta(alpha, c, t(i)))
# **************************************************
# * C a u t i o n: *
# * The theory that makes statistical *
# * diagnostics useful is derived for *
# * linear regression, with no upper- or *
# * lower-bounds on variables. *
# * The relevance of these quantities to our *
# * nonlinear model is determined by how well *
# * the linearized model (Taylor series model) *
# * eta(alpha_true, c_true) *
# * + Phi * (c - c_true) *
# * + dPhi * (alpha - alpha_true) *
# * fits the data in the neighborhood of the *
# * true values for alpha and c, where Phi *
# * and dPhi contain the partial derivatives *
# * of the model with respect to the c and *
# * alpha parameters, respectively, and are *
# * defined in ada. *
# **************************************************
#
# CorMx: (n+q) x (n+q)
# This is the estimated correlation matrix for the
# parameters. The linear parameters c are ordered
# first, followed by the nonlinear parameters alpha.
# std_dev_param: length n+q
# This vector contains the estimate of the standard
# deviation for each parameter.
# The k-th element is the square root of the k-th main
# diagonal element of the covariance matrix CovMatrix
#
#---------------------------------------------------------------
# Specification of the function ada, which computes information
# related to Phi:
#
# Phi,dPhi,Ind = ada(alpha)
#
# This function computes Phi and dPhi.
#
# On Input:
#
# alpha length q contains the current value of the alpha parameters.
#
# Note: If more input arguments are needed, call with a lambda.
# For example, if the input arguments to ada are alpha and t,
# then before calling varpro, initialize t and call varpro with
# "lambda alpha = None: ada(alpha,t)"
#
# On Output:
#
# Phi m x n1 where Phi(i,j) = phi_j(alpha,t_i).
# (n1 = n if there is no extra term;
# n1 = n+1 if an extra term is used for a nonlinear
# term with no linear coefficient)
# dPhi m x p where the columns contain partial derivative
# information for Phi and p is the number of
# columns in Ind
# Use numerical differentiation if analytical derivatives
# not available
# Ind 2 x p Column k of dPhi contains the partial
# derivative of Phi_j with respect to alpha_i,
# evaluated at the current value of alpha,
# where j = Ind(0,k) and i = Ind(1,k).
# Columns of dPhi that are always zero, independent
# of alpha, need not be stored.
# Example: If phi_0 is a function of alpha_1 and alpha_2,
# and phi_1 is a function of alpha_0 and alpha_1, then
# we can set
# Ind = [ 0 0 1 1
# 1 2 0 1 ]
# In this case, the p=4 columns of dPhi contain
# d phi_0 / d alpha_1,
# d phi_0 / d alpha_2,
# d phi_1 / d alpha_0,
# d phi_1 / d alpha_1,
# evaluated at each t_i.
# There are no restrictions on how the columns of
# dPhi are ordered, as long as Ind correctly specifies
# the ordering.
#
#---------------------------------------------------------------
#
# Note that another non-linear least squares solver can be substituted
# for the scipy one used. Also, the method "dogbox" can be replaced
# by one of the others available in the scipy routine and default
# tolerances in that routine can also be replaced.
#
# Any linear parameters that require upper or lower bounds should be put in
# alpha, not c, and treated non-linearly. (This should rarely be needed.)
#
#%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
m = len(y)
m1 = len(w)
m2 = len(t)
if (m1 != m):
raise Exception('y and w must be vectors of the same length')
if not (m == m2 or m == 2*m2):
raise Exception('y and w must be the same or twice the length of t '
'in the case where the data is complex, since in '
'this case it should be organized as (data.real, '
'data.imag)')
if (len(alpha.shape) > 1):
raise Exception('alpha must be a 1d vector containing initial guesses for nonlinear parameters')
q = len(alpha)
if (q == 0):
raise Exception('No nonlinear parameters: use linear least squares')
if (bounds is not None):
if type(bounds) is not tuple:
raise Exception('bounds must be omitted or supplied as a tuple')
q1 = len(bounds)
if (q1 != 2):
raise Exception('must specify both lower and upper bounds')
#rely on scipy least_squares for further checking of bounds
else:
bounds=(-np.inf, np.inf) #default for scipy least_squares
W = spdiags(w,0,m,m) #convert w from 1-d to 2-d array
Phi,dPhi,Ind = ada(alpha)
m1,n1 = Phi.shape
m2,n2 = dPhi.shape
ell,n3 = Ind.shape
if np.logical_and((m1 != m2),(m2 > 0)):
raise Exception('In user function ada: Phi and dPhi must have the same number of rows.')
if np.logical_or((n1 < n),(n1 > n + 1)):
raise Exception('In user function ada: The number of columns in Phi must be n or n+1.')
if np.logical_and((n2 > 0),(ell != 2)):
raise Exception('In user function ada: Ind must have two rows.')
if np.logical_and((n2 > 0),(n2 != n3)):
raise Exception('In user function ada: dPhi and Ind must have the same number of columns.')
def formJacobian(alpha, Phi, dPhi):
# indexing matches number of terms with linear coefficients
U,S,V = np.linalg.svd(W * Phi[:,:n])
if (n >= 1):
s = S #S is a vector in Python, not a matrix
else: #no linear parameters
if len(Ind)==0:
Jacobian = []
else:
Jacobian = np.zeros((len(y),len(alpha)))
Jacobian[:,Ind[2,:]] = - W * dPhi
c = []
y_est = Phi
wresid = W * (y - y_est)
myrank = 1
return Jacobian,c,wresid,y_est,myrank
#tol = m * sys.float_info.epsilon
tol = m * np.finfo(float).eps
myrank = sum(s > tol * s[0])
s = s[0:myrank]
if (myrank < n):
raise RankError(myrank,n,np.zeros(n),np.zeros(q))
yuse = y
if (n < n1):
# indexing matches number of terms with linear coefficients
yuse = y - Phi[:,-1]
temp = np.ndarray.flatten(np.transpose(U[:,0:myrank]).dot(W.dot(yuse)))
c = (temp / s).dot(V[0:myrank])
y_est = Phi[:,0:n].dot(c)
wresid = W * (yuse - y_est)
if (n < n1):
# indexing matches number of terms with linear coefficients
y_est = y_est + Phi[:,-1]
if len(dPhi)==0:
Jacobian = []
return Jacobian,c,wresid,y_est,myrank
WdPhi = W * dPhi
WdPhi_r = wresid.dot(WdPhi)
T2 = np.zeros((n1,q))
if (n1 > n): # UPDATED: should be correct
c = np.concatenate((np.array(c), np.array([1]))) # fixed concatenation here
Jac1 = np.zeros((m,q))
for j in np.arange(q):
range = np.where(Ind[1,:] == j)[0]
indrows = Ind[0,range]
Jac1[:,j] = WdPhi[:,range].dot(c[indrows])
T2[indrows,j] = WdPhi_r[range]
Jac1 = U[:,myrank:m].dot(
np.transpose(U[:,myrank:m]).dot(Jac1))
T2 = np.diag(1 / s[0:myrank]).dot(V[0:myrank]
.dot(T2[0:n,:]))
Jac2 = U[:,0:myrank].dot(T2)
Jacobian = - (Jac1 + Jac2)
return Jacobian,c,wresid,y_est,myrank # end of formJacobian
def f_lsq(alpha_trial):
Phi_trial,dPhi_trial,Ind = ada(alpha_trial)
Jacobian,c,wr_trial,y_est,myrank = formJacobian(alpha_trial,
Phi_trial,dPhi_trial)
return wr_trial,Jacobian,Phi_trial,dPhi_trial,y_est,myrank
# end of f_lsq
class Func_jacobian:
#Computes function and Jacobian in a single routine for efficiency,
#but supplies them as separate functions
def __init__(self,xold):
self.x = xold
self.fun_jac(xold)
def fun(self,x):
if not np.array_equal(self.x, x):
self.x = x
self.fun_jac(x)
return self.f
def jac(self,x):
if not np.array_equal(self.x, x):
self.x = x
self.fun_jac(x)
return self.j
def fun_jac(self,x):
wr_trial,Jacobian,Phi_trial,dPhi_trial,y_est,myrank = f_lsq(x)
self.f=wr_trial
self.j=Jacobian
fj = Func_jacobian(alpha)
result = least_squares(lambda z:fj.fun(z),alpha,lambda z: fj.jac(z),
bounds, max_nfev=max_nfev, **kwargs)
Phi,dPhi,Ind = ada(result.x)
Jacobian,c,wresid,y_est,myrank=formJacobian(result.x,Phi,dPhi)
if result.status == 0:
#maximum number of nonlinear fit iterations was exceeded
raise MaxIterationError(c, result.x)
if myrank < n:
#linear parameters are ill-determined
raise RankError(myrank, n, c, result.x)
logger.info("residual_norm={}".format(result.cost))
logger.info("gradient norm={}".format(result.optimality))
logger.info("nfev = {}".format(result.nfev))
logger.info("njev = {}".format(result.njev))
logger.info("status type = {}".format(result.status))
logger.info("status = {}".format(result.message))
wresid_norm = np.linalg.norm(wresid)
xx,pp = dPhi.shape
J = np.zeros((m,q))
for kk in np.arange(pp):
j = Ind[0,kk]
i = Ind[1,kk]
if (j > n):
J[:,i] = J[:,i] + dPhi[:,kk]
else:
J[:,i] = J[:,i] + c[j] * dPhi[:,kk]
Mat = W.dot(np.concatenate((Phi[:,np.arange(n)],J),axis=1))
Qj,Rj,Pj = linalg.qr(Mat, mode='economic', pivoting=True)
# Updated Rj matrix to have a 1 on the diagonal in the final entry corresponding to the final term
# in the model with no linear coefficient (dimension disagreement error for n != n1)
if n < n1:
np.fill_diagonal(Rj, 1)
T2 = linalg.solve_triangular(Rj,(np.identity(Rj.shape[0])))
sigma2 = wresid_norm*wresid_norm/(m-n-q)
CovMx = sigma2 * T2.dot(np.transpose(T2))
CovMatrix = np.empty(Rj.shape)
tempMat = np.empty(Rj.shape)
tempMat[:,Pj] = CovMx #Undo the pivoting permutation
CovMatrix[Pj,:] = tempMat
d = 1 / np.sqrt(np.diag(CovMatrix))
D = spdiags(d,0,n + q,n + q)
CorMx = D * CovMatrix * D
std_dev_params = np.sqrt(np.diag(CovMatrix))
logger.info('VARPRO Results:')
logger.info("linear parameters={}".format(c))
logger.info("nonlinear parameters={}".format(result.x))
logger.info('std_dev_params={}'.format(std_dev_params))
logger.info("wresid_norm={}".format(wresid_norm))
logger.info('Corr. Matrix={}\n'.format(CorMx))
return result.x,c,wresid,wresid_norm,y_est,CorMx,std_dev_params # end of varpro