2 ;;; Copyright (c) 2005--2007, by A.J. Rossini <blindglobe@gmail.com>
3 ;;; See COPYRIGHT file for any additional restrictions (BSD license).
4 ;;; Since 1991, ANSI was finally finished. Edited for ANSI Common Lisp.
6 ;;; what this should do:
7 ;;; #1 - use CFFI (and possibly Verazanno) to import C/C++.
8 ;;; #2 - what to do for Fortran? Possibly: C <-> bridge, or CLapack?
9 ;;; problem: would be better to have access to Fortran. For
10 ;;; example, think of Doug Bates comment on reverse-calls (as
11 ;;; distinct from callbacks). It would be difficult if we don't
12 ;;; -- however, has anyone run Lapack or similar through F2CL?
13 ;;; Answer: yes, Matlisp does this.
16 ;;;; linalg -- Lisp-Stat interface to basic linear algebra routines.
18 ;;;; Copyright (c) 1991, by Luke Tierney. Permission is granted for
19 ;;;; unrestricted use.
25 (defpackage :lisp-stat-linalg
29 (:export chol-decomp lu-decomp lu-solve determinant inverse sv-decomp
30 qr-decomp rcondest make-rotation spline kernel-dens kernel-smooth
31 fft make-sweep-matrix sweep-operator ax
+y numgrad numhess
35 (in-package #:lisp-stat-linalg
)
38 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
40 ;;;; Lisp to C number conversion and checking
42 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
45 ;;;; Lisp to/from C sequence and matrix conversion and checking
49 "FIXME:AJR this is not used anywhere?"
52 (defun check-fixnum (a)
53 (if (/= 0 (la-data-mode a
)) (error "not an integer sequence - ~s" a
)))
55 (defun check-real (data)
56 (let ((data (compound-data-seq data
)))
59 (let ((n (length data
)))
63 (check-one-real (aref data i
)))))
64 ((consp data
) (dolist (x data
) (check-one-real x
)))
65 (t (error "bad sequence - ~s" data
)))))
67 (defun vec-assign (a i x
) (setf (aref a i
) x
))
69 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
70 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
72 ;;;; Lisp Interfaces to Linear Algebra Routines
74 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
75 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
77 ;;; FIXME: use dpbt[f2|rf], dpbstf, dpot[f2|rf]; dpptrf, zpbstf, zpbt[f2|rf]
78 ;;; remember: factorization = decomposition, depending on training.
80 (defun chol-decomp (a &optional
(maxoffl 0.0))
82 Modified Cholesky decomposition. A should be a square, symmetric matrix.
83 Computes lower triangular matrix L such that L L^T = A + D where D is a
84 diagonal matrix. If A is strictly positive definite D will be zero.
85 Otherwise, D is as small as possible to make A + D numerically strictly
86 positive definite. Returns a list (L (max D))."
87 (check-square-matrix a
)
89 (let* ((n (array-dimension a
0))
90 (result (make-array (list n n
)))
91 (dpars (list maxoffl
0.0)))
93 (let ((mat (la-data-to-matrix a
+mode-re
+))
94 (dp (la-data-to-vector dpars
+mode-re
+)))
97 (chol-decomp-front mat n dp
)
98 (la-matrix-to-data mat n n
+mode-re
+ result
)
99 (la-vector-to-data dp
2 +mode-re
+ dpars
))
100 (la-free-matrix mat n
)
101 (la-free-vector dp
)))
102 (list result
(second dpars
))))
107 ;;; i.e. result use by:
108 ;;; (setf (values (lu-out1 lu-out2 lu-out3)) (matlisp:lu my-matrix))
109 ;;; for solution, ...
111 ;;; (matlisp:gesv a b &opt ipivot)
115 A is a square matrix of numbers (real or complex). Computes the LU
116 decomposition of A and returns a list of the form (LU IV D FLAG), where
117 LU is a matrix with the L part in the lower triangle, the U part in the
118 upper triangle (the diagonal entries of L are taken to be 1), IV is a vector
119 describing the row permutation used, D is 1 if the number of permutations
120 is odd, -1 if even, and FLAG is T if A is numerically singular, NIL otherwise.
122 (check-square-matrix a
)
123 (let* ((n (array-dimension a
0))
124 (mode (max +mode-re
+ (la-data-mode a
)))
125 (result (list (make-array (list n n
)) (make-array n
) nil nil
)))
126 (let ((mat (la-data-to-matrix a mode
))
127 (iv (la-vector n
+mode-in
+))
128 (d (la-vector 1 +mode-re
+))
132 (setf singular
(lu-decomp-front mat n iv mode d
))
133 (la-matrix-to-data mat n n mode
(first result
))
134 (la-vector-to-data iv n
+mode-in
+ (second result
))
135 (setf (third result
) (la-get-double d
0))
136 (setf (fourth result
) (if (= singular
0.0) nil t
)))
137 (la-free-matrix mat n
)
142 (defun lu-solve (lu lb
)
144 LU is the result of (LU-DECOMP A) for a square matrix A, B is a sequence.
145 Returns the solution to the equation Ax = B. Signals an error if A is
147 (let ((la (first lu
))
149 (check-square-matrix la
)
150 (check-sequence lidx
)
153 (let* ((n (num-rows la
))
154 (result (make-sequence (if (consp lb
) 'list
'vector
) n
))
155 (a-mode (la-data-mode la
))
156 (b-mode (la-data-mode lb
)))
157 (if (/= n
(length lidx
)) (error "index sequence is wrong length"))
158 (if (/= n
(length lb
)) (error "right hand side is wrong length"))
159 (let* ((mode (max +mode-re
+ a-mode b-mode
))
160 (a (la-data-to-matrix la mode
))
161 (indx (la-data-to-vector lidx
+mode-in
+))
162 (b (la-data-to-vector lb mode
))
166 (setf singular
(lu-solve-front a n indx b mode
))
167 (la-vector-to-data b n mode result
))
169 (la-free-vector indx
)
171 (if (/= 0.0 singular
) (error "matrix is (numerically) singular"))
174 (defun determinant (a)
176 Returns the determinant of the square matrix M."
177 (let* ((lu (lu-decomp a
))
183 (flet ((fabs (x) (float (abs x
) 0.d0
)))
184 (dotimes (i n
(* d1
(exp d2
)))
186 (let* ((x (aref la i i
))
188 (if (= 0.0 magn
) (return 0.d0
))
189 (setf d1
(* d1
(/ x magn
)))
190 (setf d2
(+ d2
(log magn
))))))))
194 Returns the inverse of the the square matrix M; signals an error if M is ill
195 conditioned or singular"
196 (check-square-matrix a
)
197 (let ((n (num-rows a
))
198 (mode (max +mode-re
+ (la-data-mode a
))))
200 (let ((result (make-array (list n n
) :initial-element
0)))
203 (setf (aref result i i
) 1))
204 (let ((mat (la-data-to-matrix a mode
))
205 (inv (la-data-to-matrix result mode
))
206 (iv (la-vector n
+mode-in
+))
207 (v (la-vector n mode
))
211 (setf singular
(lu-inverse-front mat n iv v mode inv
))
212 (la-matrix-to-data inv n n mode result
))
213 (la-free-matrix mat n
)
214 (la-free-matrix inv n
)
217 (if (/= singular
0) (error "matrix is (numerically) singular"))
221 ;;;; SV Decomposition
226 A is a matrix of real numbers with at least as many rows as columns.
227 Computes the singular value decomposition of A and returns a list of the form
228 (U W V FLAG) where U and V are matrices whose columns are the left and right
229 singular vectors of A and W is the sequence of singular values of A. FLAG is T
230 if the algorithm converged, NIL otherwise."
232 (let* ((m (num-rows a
))
234 (mode (max +mode-re
+ (la-data-mode a
)))
235 (result (list (make-array (list m n
))
237 (make-array (list n n
))
239 (if (< m n
) (error "number of rows less than number of columns"))
240 (if (= mode
+mode-cx
+) (error "complex SVD not available yet"))
241 (let ((mat (la-data-to-matrix a mode
))
242 (w (la-vector n
+mode-re
+))
243 (v (la-matrix n n
+mode-re
+))
247 (setf converged
(sv-decomp-front mat m n w v
))
248 (la-matrix-to-data mat m n mode
(first result
))
249 (la-vector-to-data w n mode
(second result
))
250 (la-matrix-to-data v n n mode
(third result
))
251 (setf (fourth result
) (if (/= 0.0 converged
) t nil
)))
252 (la-free-matrix mat m
)
254 (la-free-matrix v n
))
259 ;;;; QR Decomposition
262 (defun qr-decomp (a &optional pivot
)
263 "Args: (a &optional pivot)
264 A is a matrix of real numbers with at least as many rows as columns. Computes
265 the QR factorization of A and returns the result in a list of the form (Q R).
266 If PIVOT is true the columns of X are first permuted to place insure the
267 absolute values of the diagonal elements of R are nonincreasing. In this case
268 the result includes a third element, a list of the indices of the columns in
269 the order in which they were used."
271 (let* ((m (num-rows a
))
273 (mode (max +mode-re
+ (la-data-mode a
)))
276 (list (make-array (list m n
))
277 (make-array (list n n
))
279 (list (make-array (list m n
)) (make-array (list n n
))))))
280 (if (< m n
) (error "number of rows less than number of columns"))
281 (if (= mode
+mode-cx
+) (error "complex QR decomposition not available yet"))
282 (let ((mat (la-data-to-matrix a mode
))
283 (v (la-matrix n n
+mode-re
+))
284 (jpvt (la-vector n
+mode-in
+)))
287 (qr-decomp-front mat m n v jpvt p
)
288 (la-matrix-to-data mat m n mode
(first result
))
289 (la-matrix-to-data v n n mode
(second result
))
290 (if pivot
(la-vector-to-data jpvt n
+mode-in
+ (third result
))))
291 (la-free-matrix mat m
)
293 (la-free-vector jpvt
))
297 ;;;; Estimate of Condition Number for Lower Triangular Matrix
302 Returns an estimate of the reciprocal of the L1 condition number of an upper
303 triangular matrix a."
304 (check-square-matrix a
)
305 (let ((mode (max +mode-re
+ (la-data-mode a
)))
307 (if (= mode
+mode-cx
+)
308 (error "complex condition estimate not available yet"))
309 (let ((mat (la-data-to-matrix a mode
))
312 (setf est
(rcondest-front mat n
))
313 (la-free-matrix mat n
))
317 ;;;; Make Rotation Matrix
320 (defun make-rotation (x y
&optional alpha
)
321 "Args: (x y &optional alpha)
322 Returns a rotation matrix for rotating from X to Y, or from X toward Y
323 by angle ALPHA, in radians. X and Y are sequences of the same length."
326 (if alpha
(check-one-real alpha
))
327 (let* ((n (length x
))
328 (mode (max +mode-re
+ (la-data-mode x
) (la-data-mode y
)))
329 (use-angle (if alpha
1 0))
330 (angle (if alpha
(float alpha
0.0) 0.0))
331 (result (make-array (list n n
))))
332 (if (/= n
(length y
)) (error "sequences not the same length"))
333 (if (= mode
+mode-cx
+) (error "complex data not supported yet"))
334 (let ((px (la-data-to-vector x
+mode-re
+))
335 (py (la-data-to-vector y
+mode-re
+))
336 (rot (la-matrix n n
+mode-re
+)))
339 (make-rotation-front n rot px py use-angle angle
)
340 (la-matrix-to-data rot n n
+mode-re
+ result
))
343 (la-free-matrix rot n
))
347 ;;;; Eigenvalues and Vectors
352 Returns list of list of eigenvalues and list of eigenvectors of square,
353 symmetric matrix A. Third element of result is NIL if algorithm converges.
354 If the algorithm does not converge, the third element is an integer I.
355 In this case the eigenvalues 0, ..., I are not reliable."
356 (check-square-matrix a
)
357 (let ((mode (max +mode-re
+ (la-data-mode a
)))
359 (if (= mode
+mode-cx
+) (error "matrix must be real and symmetric"))
360 (let ((evals (make-array n
))
361 (evecs (make-list (* n n
)))
362 (pa (la-data-to-vector (compound-data-seq a
) +mode-re
+))
363 (w (la-vector n
+mode-re
+))
364 (z (la-vector (* n n
) +mode-re
+))
365 (fv1 (la-vector n
+mode-re
+))
369 (setf ierr
(eigen-front pa n w z fv1
))
370 (la-vector-to-data w n
+mode-re
+ evals
)
371 (la-vector-to-data z
(* n n
) +mode-re
+ evecs
))
375 (la-free-vector fv1
))
376 (list (nreverse evals
)
377 (nreverse (mapcar #'(lambda (x) (coerce x
'vector
))
378 (split-list evecs n
)))
379 (if (/= 0 ierr
) (- n ierr
))))))
382 ;;;; Spline Interpolation
385 (defun make-smoother-args (x y xvals
)
391 (unless (integerp xvals
)
392 (check-sequence xvals
)
394 (let* ((n (length x
))
395 (ns (if (integerp xvals
) xvals
(length xvals
)))
396 (result (list (make-list ns
) (make-list ns
))))
397 (if (and y
(/= n
(length y
))) (error "sequences not the same length"))
398 (list x y n
(if (integerp xvals
) 0 1) ns xvals result
)))
400 (defun get-smoother-result (args) (seventh args
))
402 (defmacro with-smoother-data
((x y xvals is-reg
) &rest body
)
409 (unless (integerp ,xvals
)
410 (check-sequence ,xvals
)
412 (let* ((supplied (not (integerp ,xvals
)))
414 (ns (if supplied
(length ,xvals
) ,xvals
))
415 (result (list (make-list ns
) (make-list ns
))))
416 (if (and ,is-reg
(/= n
(length ,y
)))
417 (error "sequences not the same length"))
418 (if (and (not supplied
) (< ns
2))
419 (error "too few points for interpolation"))
420 (let* ((px (la-data-to-vector ,x
+mode-re
+))
421 (py (if ,is-reg
(la-data-to-vector ,y
+mode-re
+)))
423 (la-data-to-vector ,xvals
+mode-re
+)
424 (la-vector ns
+mode-re
+)))
425 (pys (la-vector ns
+mode-re
+)))
426 (unless supplied
(la-range-to-rseq n px ns pxs
))
429 (la-vector-to-data pxs ns
+mode-re
+ (first result
))
430 (la-vector-to-data pys ns
+mode-re
+ (second result
)))
432 (if ,is-reg
(la-free-vector py
))
434 (la-free-vector pys
))
437 (defun spline (x y
&key
(xvals 30))
438 "Args: (x y &key xvals)
439 Returns list of x and y values of natural cubic spline interpolation of (X,Y).
440 X must be strictly increasing. XVALS can be an integer, the number of equally
441 spaced points to use in the range of X, or it can be a sequence of points at
442 which to interpolate."
443 (with-smoother-data (x y xvals t
)
444 (let ((work (la-vector (* 2 n
) +mode-re
+))
447 (setf error
(spline-front n px py ns pxs pys work
))
448 (la-free-vector work
))
449 (if (/= error
0) (error "bad data for splines")))))
452 ;;;; Kernel Density Estimators and Smoothers
455 (defun kernel-type-code (type)
456 (cond ((eq type
'u
) 0)
461 (defun kernel-dens (x &key
(type 'b
) (width -
1.0) (xvals 30))
462 "Args: (x &key xvals width type)
463 Returns list of x and y values of kernel density estimate of X. XVALS can be an
464 integer, the number of equally spaced points to use in the range of X, or it
465 can be a sequence of points at which to interpolate. WIDTH specifies the
466 window width. TYPE specifies the lernel and should be one of the symbols G, T,
467 U or B for gaussian, triangular, uniform or bisquare. The default is B."
468 (check-one-real width
)
469 (with-smoother-data (x nil xvals nil
) ;; warning about deleting unreachable code is TRUE -- 2nd arg=nil!
470 (let ((code (kernel-type-code type
))
472 (setf error
(kernel-dens-front px n width pxs pys ns code
))
473 (if (/= 0 error
) (error "bad kernel density data")))))
475 (defun kernel-smooth (x y
&key
(type 'b
) (width -
1.0) (xvals 30))
476 "Args: (x y &key xvals width type)
477 Returns list of x and y values of kernel smooth of (X,Y). XVALS can be an
478 integer, the number of equally spaced points to use in the range of X, or it
479 can be a sequence of points at which to interpolate. WIDTH specifies the
480 window width. TYPE specifies the lernel and should be one of the symbols G, T,
481 U or B for Gaussian, triangular, uniform or bisquare. The default is B."
482 (check-one-real width
)
483 (with-smoother-data (x y xvals t
)
484 (let ((code (kernel-type-code type
))
486 (kernel-smooth-front px py n width pxs pys ns code
)
487 ;; if we get the Lisp version ported from C, uncomment below and
488 ;; comment above. (thanks to Carlos Ungil for the initial CFFI
490 ;;(kernel-smooth-Cport px py n width pxs pys ns code)
491 (if (/= 0 error
) (error "bad kernel density data")))))
495 (defun kernel-smooth-Cport (px py n width
;;wts wds ;; see above for mismatch?
497 "Port of kernel_smooth (Lib/kernel.c) to Lisp.
498 FIXME:kernel-smooth-Cport : This is broken.
499 Until this is fixed, we are using Luke's C code and CFFI as glue."
501 ((and (< n
2) (<= width
0)) 1.0)
502 (t (let* ((xmin (min px
))
504 (width (/ (- xmax xmin
) (+ 1.0 (log n
)))))
505 (dotimes (i (- ns
1))
509 (dotimes (j (- n
1)) )
510 ;;;possible nasty errors...
513 ((lwidth (if wds
(* width
(aref wds j
)) width
))
514 (lwt (* (kernel-Cport (aref xs i
) (aref px j
) lwidth ktype
) ;; px?
515 (if wts
(aref wts j
) 1.0))))
516 (setf wsum
(+ wsum lwt
))
517 (setf ysum
(if py
(+ ysum
(* lwt
(aref py j
)))))) ;; py? y?
527 (defun kernel-Cport (x y w ktype
)
528 "Port of kernel() (Lib/kernel.c) to Lisp.
529 x,y,w are doubles, type is an integer"
533 (cond ((eq ktype
"B")
538 (/ (/ (* 15.0 (* (- 1.0 (* 4 z z
)) ;; k/w
539 (- 1.0 (* 4 z z
)))) ;; k/w
544 (let* ((w (* w
0.25))
546 (k (/ (exp (* -
0.5 z z
))
552 (k (if (< (abs z
) 0.5)
557 (cond ((and (> z -
1.0)
568 ;;;; Lowess Smoother Interface
571 (defun |base-lowess|
(s1 s2 f nsteps delta
)
577 (check-one-fixnum nsteps
)
578 (check-one-real delta
)
579 (let* ((n (length s1
))
580 (result (make-list n
)))
581 (if (/= n
(length s2
)) (error "sequences not the same length"))
582 (let ((x (la-data-to-vector s1
+mode-re
+))
583 (y (la-data-to-vector s2
+mode-re
+))
584 (ys (la-vector n
+mode-re
+))
585 (rw (la-vector n
+mode-re
+))
586 (res (la-vector n
+mode-re
+))
590 (setf error
(base-lowess-front x y n f nsteps delta ys rw res
))
591 (la-vector-to-data ys n
+mode-re
+ result
))
596 (la-free-vector res
))
597 (if (/= error
0) (error "bad data for lowess"))
601 static LVAL add_contour_point
(i, j
, k
, l
, x
, y
, z
, v
, result
)
611 if
((z[i][j] <= v && v < z[k][l]) || (z[k][l] <= v && v < z[i][j])) {
614 p = (v - z[i][j]) / (z[k][l] - z[i][j]);
616 rplaca(pt, cvflonum((FLOTYPE) (q * x[i] + p * x[k])));
617 rplaca
(cdr(pt), cvflonum
((FLOTYPE) (q * y
[j] + p * y[l])));
618 result = cons(pt, result);
624 LVAL xssurface_contour()
626 LVAL s1, s2, mat, result;
632 s1 = xsgetsequence();
633 s2 = xsgetsequence();
635 v = makedouble(xlgetarg());
638 n = seqlen(s1); m = seqlen(s2);
639 if (n != numrows(mat) || m != numcols(mat)) xlfail("dimensions do not match");
640 if (data_mode(s1) == CX || data_mode(s2) == CX || data_mode(mat) == CX)
641 xlfail("data must be real");
643 x = (RVector) data_to_vector(s1, RE);
644 y = (RVector) data_to_vector(s2, RE);
645 z = (RMatrix) data_to_matrix(mat, RE);
649 for (i = 0; i < n - 1; i++) {
650 for (j = 0; j < m - 1; j++) {
651 result = add_contour_point(i, j, i, j+1, x, y, z, v, result);
652 result = add_contour_point(i, j+1, i+1, j+1, x, y, z, v, result);
653 result = add_contour_point(i+1, j+1, i+1, j, x, y, z, v, result);
654 result = add_contour_point(i+1, j, i, j, x, y, z, v, result);
671 ;;; ??replace with matlisp:fft and matlisp:ifft (the latter for inverse mapping)
673 (defun fft (x &optional inverse)
674 "Args: (x &optional inverse)
675 Returns unnormalized Fourier transform of X, or inverse transform if INVERSE
678 (let* ((n (length x))
679 (mode (la-data-mode x))
680 (isign (if inverse -1 1))
681 (result (if (consp x) (make-list n) (make-array n))))
682 (let ((px (la-data-to-vector x +mode-cx+))
683 (work (la-vector (+ (* 4 n) 15) +mode-re+)))
686 (fft-front n px work isign)
687 (la-vector-to-data px n +mode-cx+ result))
689 (la-free-vector work))
693 ;;; SWEEP Operator: FIXME: use matlisp
696 (defun make-sweep-front (x y w n p mode has_w x_mean result)
697 (declare (fixnum n p mode has_w))
712 (has-w (if (/= 0 has_w) t nil))
714 (declare (long-float val dxi dyi dv dw sum_w dxik dxjk dyj
715 dx_meani dx_meanj dy_mean))
716 ;; (declare-double val dxi dyi dv dw sum_w dxik dxjk dyj
717 ;; dx_meani dx_meanj dy_mean)
719 (if (> mode RE) (error "not supported for complex data yet"))
721 (setf x_data (compound-data-seq x))
722 (setf result_data (compound-data-seq result))
724 ;; find the mean of y
729 (setf dyi (makedouble (aref y i)))
731 (setf dw (makedouble (aref w i)))
733 (setf dyi (* dyi dw)))
735 (if (not has-w) (setf sum_w (float n 0.0)))
736 (if (<= sum_w 0.0) (error "non positive sum of weights"))
737 (setf dy_mean (/ val sum_w))
739 ;; find the column means
745 (setf dxi (makedouble (aref x_data (+ (* p i) j))))
747 (setf dw (makedouble (aref w i)))
748 (setf dxi (* dxi dw)))
750 (setf (aref x_mean j) (/ val sum_w)))
752 ;; put 1/sum_w in topleft, means on left, minus means on top
753 (setf (aref result_data 0) (/ 1.0 sum_w))
756 (setf dxi (makedouble (aref x_mean i)))
757 (setf (aref result_data (+ i 1)) (- dxi))
758 (setf (aref result_data (* (+ i 1) (+ p 2))) dxi))
759 (setf (aref result_data (+ p 1)) (- dy_mean))
760 (setf (aref result_data (* (+ p 1) (+ p 2))) dy_mean)
762 ;; put sums of adjusted cross products in body
770 (setf dxik (makedouble (aref x_data (+ (* p k) i))))
771 (setf dxjk (makedouble (aref x_data (+ (* p k) j))))
772 (setf dx_meani (makedouble (aref x_mean i)))
773 (setf dx_meanj (makedouble (aref x_mean j)))
774 (setf dv (* (- dxik dx_meani) (- dxjk dx_meanj)))
776 (setf dw (makedouble (aref w k)))
779 (setf (aref result_data (+ (* (+ i 1) (+ p 2)) (+ j 1))) val)
780 (setf (aref result_data (+ (* (+ j 1) (+ p 2)) (+ i 1))) val))
784 (setf dxik (makedouble (aref x_data (+ (* p j) i))))
785 (setf dyj (makedouble (aref y j)))
786 (setf dx_meani (makedouble (aref x_mean i)))
787 (setf dv (* (- dxik dx_meani) (- dyj dy_mean)))
789 (setf dw (makedouble (aref w j)))
792 (setf (aref result_data (+ (* (+ i 1) (+ p 2)) (+ p 1))) val)
793 (setf (aref result_data (+ (* (+ p 1) (+ p 2)) (+ i 1))) val))
797 (setf dyj (makedouble (aref y j)))
798 (setf dv (* (- dyj dy_mean) (- dyj dy_mean)))
800 (setf dw (makedouble (aref w j)))
803 (setf (aref result_data (+ (* (+ p 1) (+ p 2)) (+ p 1))) val)))
805 ;;; FIXME: use matlisp
806 (defun sweep-in-place-front (a rows cols mode k tol)
807 "Sweep algorithm for linear regression."
808 (declare (long-float tol))
809 (declare (fixnum rows cols mode k))
817 (declare (long-float pivot aij aik akj akk))
819 (if (> mode RE) (error "not supported for complex data yet"))
820 (if (or (< k 0) (>= k rows) (>= k cols)) (error "index out of range"))
822 (setf tol (max tol machine-epsilon))
823 (setf data (compound-data-seq a))
825 (setf pivot (makedouble (aref data (+ (* cols k) k))))
828 ((or (> pivot tol) (< pivot (- tol)))
833 (when (and (/= i k) (/= j k))
834 (setf aij (makedouble (aref data (+ (* cols i) j))))
835 (setf aik (makedouble (aref data (+ (* cols i) k))))
836 (setf akj (makedouble (aref data (+ (* cols k) j))))
837 (setf aij (- aij (/ (* aik akj) pivot)))
838 (setf (aref data (+ (* cols i) j)) aij))))
842 (setf aik (makedouble (aref data (+ (* cols i) k))))
844 (setf aik (/ aik pivot))
845 (setf (aref data (+ (* cols i) k)) aik)))
849 (setf akj (makedouble (aref data (+ (* cols k) j))))
851 (setf akj (- (/ akj pivot)))
852 (setf (aref data (+ (* cols k) j)) akj)))
854 (setf akk (/ 1.0 pivot))
855 (setf (aref data (+ (* cols k) k)) akk)
859 ;; FIXME: use matlisp
860 (defun make-sweep-matrix (x y &optional w)
861 "Args: (x y &optional weights)
862 X is matrix, Y and WEIGHTS are sequences. Returns the sweep matrix of the
863 (weighted) regression of Y on X"
866 (if w (check-sequence w))
867 (let ((n (num-rows x))
869 (if (/= n (length y)) (error "dimensions do not match"))
870 (if (and w (/= n (length w))) (error "dimensions do not match"))
871 (let ((mode (max (la-data-mode x)
873 (if w (la-data-mode w) 0)))
874 (result (make-array (list (+ p 2) (+ p 2))))
875 (x-mean (make-array p))
876 (y (coerce y 'vector))
877 (w (if w (coerce w 'vector)))
879 (make-sweep-front x y w n p mode has-w x-mean result)
882 (defun sweep-in-place (a k tol)
886 (let ((rows (num-rows a))
888 (mode (la-data-mode a)))
889 (let ((swept (sweep-in-place-front a rows cols mode k tol)))
890 (if (/= 0 swept) t nil))))
892 (defun sweep-operator (a columns &optional tolerances)
893 "Args: (a indices &optional tolerances)
894 A is a matrix, INDICES a sequence of the column indices to be swept. Returns
895 a list of the swept result and the list of the columns actually swept. (See
896 MULTREG documentation.) If supplied, TOLERANCES should be a list of real
897 numbers the same length as INDICES. An index will only be swept if its pivot
898 element is larger than the corresponding element of TOLERANCES."
900 (check-sequence columns)
901 (if tolerances (check-sequence tolerances))
903 (check-fixnum columns)
904 (if tolerances (check-real tolerances))
906 (result (copy-array a))
908 (columns (coerce columns 'list) (cdr columns))
909 (tolerances (if (consp tolerances) (coerce tolerances 'list))
910 (if (consp tolerances) (cdr tolerances))))
911 ((null columns) (list result swept-columns))
912 (let ((col (first columns))
913 (tol (if (consp tolerances) (first tolerances) tol)))
914 (if (sweep-in-place result col tol)
915 (setf swept-columns (cons col swept-columns))))))
924 (defun ax+y (a x y &optional lower)
925 "Args (a x y &optional lower)
926 Returns (+ (matmult A X) Y). If LOWER is not nil, A is taken to be lower
928 This could probably be made more efficient."
929 (check-square-matrix a)
935 (let* ((n (num-rows a))
936 (result (make-list n))
937 (a (compound-data-seq a)))
939 (if (or (/= n (length x)) (/= n (length y)))
940 (error "dimensions do not match"))
941 (do* ((tx (make-next-element x) (make-next-element x))
942 (ty (make-next-element y))
943 (tr (make-next-element result))
945 (start 0 (+ start n))
946 (end (if lower (+ i 1) n) (if lower (+ i 1) n)))
948 (declare (fixnum i start end))
949 (let ((val (get-next-element ty i)))
952 (setf val (+ val (* (get-next-element tx j)
953 (aref a (+ start j))))))
954 (set-next-element tr i val)))))
957 ;;;; Maximization and Numerical Derivatives
960 (defvar *maximize-callback-function* nil)
961 (defvar *maximize-callback-arg* nil)
963 (defun data2double (n data ptr)
965 (let* ((seq (compound-data-seq data))
966 (elem (make-next-element seq)))
967 (if (/= (length seq) n) (error "bad data size"))
970 (la-put-double ptr i (get-next-element elem i)))))
972 (defun maximize-callback (n px pfval pgrad phess pderivs)
973 (la-vector-to-data px n +mode-re+ *maximize-callback-arg*)
974 (let* ((val (funcall *maximize-callback-function* *maximize-callback-arg*))
975 (derivs (if (consp val) (- (length val) 1) 0)))
976 (la-put-integer pderivs 0 derivs)
977 (la-put-double pfval 0 (if (consp val) (first val) val))
978 (if (<= 1 derivs) (data2double n (second val) pgrad))
979 (if (<= 2 derivs) (data2double (* n n) (third val) phess))))
981 (defun numgrad (f x &optional scale (h -1.0))
982 "Args: (f x &optional scale derivstep)
983 Computes the numerical gradient of F at X."
987 (check-sequence scale)
990 (let* ((n (length x))
991 (result (make-list n)))
992 (if (and scale (/= n (length scale)))
993 (error "scale not the same length as x"))
994 (let ((*maximize-callback-function* f)
995 (*maximize-callback-arg* (make-list n)))
996 (let ((px (la-data-to-vector x +mode-re+))
997 (pgrad (la-vector n +mode-re+))
998 (pscale (la-data-to-vector
999 (if scale scale (make-list n :initial-element 1.0))
1003 (numgrad-front n px pgrad h pscale)
1004 (la-vector-to-data pgrad n +mode-re+ result))
1006 (la-free-vector pgrad)
1007 (la-free-vector pscale))))
1010 (defun numhess (f x &optional scale (h -1.0) all)
1011 "Args: (f x &optional scale derivstep)
1012 Computes the numerical Hessian matrix of F at X."
1016 (check-sequence scale)
1019 (let* ((n (length x))
1021 (list nil (make-list n) (make-array (list n n)))
1022 (make-array (list n n)))))
1023 (if (and scale (/= n (length scale)))
1024 (error "scale not the same length as x"))
1025 (let ((*maximize-callback-function* f)
1026 (*maximize-callback-arg* (make-list n)))
1027 (let ((hess-data (compound-data-seq (if all (third result) result)))
1028 (px (la-data-to-vector x +mode-re+))
1029 (pf (la-vector 1 +mode-re+))
1030 (pgrad (la-vector n +mode-re+))
1031 (phess (la-vector (* n n) +mode-re+))
1032 (pscale (la-data-to-vector
1033 (if scale scale (make-list n :initial-element 1.0))
1037 (numhess-front n px pf pgrad phess h pscale)
1039 (setf (first result) (la-get-double pf 0))
1040 (la-vector-to-data pgrad n +mode-re+ (second result)))
1041 (la-vector-to-data phess (* n n) +mode-re+ hess-data))
1044 (la-free-vector pgrad)
1045 (la-free-vector phess)
1046 (la-free-vector pscale))))
1049 (defun init-minfo-ipar-values (n ipars)
1061 (setf (aref ipars 0) n)
1062 (setf (aref ipars 1) m)
1063 (setf (aref ipars 2) k)
1064 (setf (aref ipars 3) itnlimit)
1065 (setf (aref ipars 4) backtrack)
1066 (setf (aref ipars 5) verbose)
1067 (setf (aref ipars 6) vals_suppl)
1068 (setf (aref ipars 7) exptilt)
1069 (setf (aref ipars 8) count)
1070 (setf (aref ipars 9) termcode)))
1072 (defun init-minfo-dpar-values (h dpars)
1081 (setf (aref dpars 0) typf)
1082 (setf (aref dpars 1) h)
1083 (setf (aref dpars 2) gradtol)
1084 (setf (aref dpars 3) steptol)
1085 (setf (aref dpars 4) maxstep)
1086 (setf (aref dpars 5) dflt)
1087 (setf (aref dpars 6) tilt)
1088 (setf (aref dpars 7) newtilt)
1089 (setf (aref dpars 8) hessadd)))
1091 (defun init-minfo-internals (n h internals)
1092 (let ((ipars (aref internals 8))
1093 (dpars (aref internals 9)))
1094 (init-minfo-ipar-values n ipars)
1095 (init-minfo-dpar-values h dpars)))
1097 (defun new-minfo-internals (f x &key scale ((:derivstep h) -1.0))
1101 (let ((n (length x)))
1103 (check-sequence scale)
1105 (if (/= n (length scale)) (error "scale and x not the same length")))
1106 (let ((internals (make-array 12)))
1107 (setf (aref internals 0) f)
1108 (setf (aref internals 3) (if (consp x) (copy-list x) (coerce x 'list)))
1109 (setf (aref internals 4)
1110 (if scale (copy-seq scale) (make-array n :initial-element 1.0)))
1111 (setf (aref internals 5) (make-list (+ 1 n (* n n))))
1112 (setf (aref internals 8) (make-array 10))
1113 (setf (aref internals 9) (make-array 9))
1114 (init-minfo-internals n h internals)
1117 (defun minfo-maximize (internals &optional verbose)
1118 "This function does what?"
1119 (let* ((f (aref internals 0))
1120 (x (aref internals 3))
1121 (fvals (aref internals 5))
1123 (v (if verbose (if (integerp verbose) verbose 1) -1)))
1124 (setf (aref internals 3) (copy-list x))
1125 (setf (aref internals 5) (copy-list fvals))
1126 (let ((*maximize-callback-function* f)
1127 (*maximize-callback-arg* (make-list n)))
1128 (let* ((x (aref internals 3))
1129 (scale (aref internals 4))
1130 (fvals (aref internals 5))
1131 (ip (aref internals 8))
1132 (dp (aref internals 9))
1133 (px (la-data-to-vector x +mode-re+))
1134 (pscale (la-data-to-vector scale +mode-re+))
1135 (pfvals (la-vector (length fvals) +mode-re+))
1136 (pip (la-data-to-vector ip +mode-in+))
1137 (pdp (la-data-to-vector dp +mode-re+)))
1140 (base-minfo-maximize px pfvals pscale pip pdp v)) ;; access to C
1141 (la-vector-to-data px n +mode-re+ x)
1142 (la-vector-to-data pfvals (+ 1 n (* n n)) +mode-re+ fvals)
1143 (la-vector-to-data pip (length ip) +mode-in+ ip)
1144 (la-vector-to-data pdp (length dp) +mode-re+ dp))
1148 ;;;; Miscellaneous Routines
1152 (defun split-list (x n)
1154 Returns a list of COLS lists of equal length of the elements of LIST.
1155 Example: (split-list '(1 2 3 4 5 6) 2) returns ((1 2 3) (4 5 6))"
1156 (check-one-fixnum n)
1157 (if (/= (rem (length x) n) 0) (error "length not divisible by ~a" n))
1158 (flet ((next-split ()
1161 (dotimes (i n result)
1162 (declare (fixnum i))
1163 (let ((c-elem (list (first x))))
1164 (cond ((null result)
1165 (setf result c-elem)
1168 (setf (rest end) c-elem)
1169 (setf end (rest end)))))
1170 (setf x (rest x))))))
1173 (k (/ (length x) n)))
1174 (declare (fixnum k))
1175 (dotimes (i k result)
1176 (declare (fixnum i))
1177 (let ((c-sub (list (next-split))))
1178 (cond ((null result)
1182 (setf (rest end) c-sub)
1183 (setf end (rest end)))))))))