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.
15 ;;; #3 - Use a lisp-based matrix system drop-in? (matlisp, femlisp, clem, ...?)
19 ;;;; linalg -- Lisp-Stat interface to basic linear algebra routines.
21 ;;;; Copyright (c) 1991, by Luke Tierney. Permission is granted for
22 ;;;; unrestricted use.
28 (defpackage :lisp-stat-linalg
34 (:shadowing-import-from
:lisp-stat-math
35 expt
+ -
* / ** mod rem abs
1+ 1- log exp sqrt sin cos tan
36 asin acos atan sinh cosh tanh asinh acosh atanh float random
37 truncate floor ceiling round minusp zerop plusp evenp oddp
38 < <= = /= >= > complex conjugate realpart imagpart phase
39 min max logand logior logxor lognot ffloor fceiling
40 ftruncate fround signum cis
)
41 (:export chol-decomp lu-decomp lu-solve determinant inverse sv-decomp
42 qr-decomp rcondest make-rotation spline kernel-dens kernel-smooth
43 fft make-sweep-matrix sweep-operator ax
+y eigen
))
45 (in-package #:lisp-stat-linalg
)
47 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
49 ;;;; Lisp to C number conversion and checking
51 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
54 ;;;; Lisp to/from C sequence and matrix conversion and checking
58 "FIXME:AJR this is not used anywhere?"
61 (defun check-fixnum (a)
62 (if (/= 0 (la-data-mode a
)) (error "not an integer sequence - ~s" a
)))
64 (defun check-real (data)
65 (let ((data (compound-data-seq data
)))
68 (let ((n (length data
)))
72 (check-one-real (aref data i
)))))
73 ((consp data
) (dolist (x data
) (check-one-real x
)))
74 (t (error "bad sequence - ~s" data
)))))
76 (defun vec-assign (a i x
) (setf (aref a i
) x
))
78 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
79 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
81 ;;;; Lisp Interfaces to Linear Algebra Routines
83 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
84 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
86 ;;; FIXME: use dpbt[f2|rf], dpbstf, dpot[f2|rf]; dpptrf, zpbstf, zpbt[f2|rf]
87 ;;; remember: factorization = decomposition, depending on training.
89 (defun chol-decomp (a &optional
(maxoffl 0.0))
91 Modified Cholesky decomposition. A should be a square, symmetric matrix.
92 Computes lower triangular matrix L such that L L^T = A + D where D is a
93 diagonal matrix. If A is strictly positive definite D will be zero.
94 Otherwise, D is as small as possible to make A + D numerically strictly
95 positive definite. Returns a list (L (max D))."
96 (check-square-matrix a
)
98 (let* ((n (array-dimension a
0))
99 (result (make-array (list n n
)))
100 (dpars (list maxoffl
0.0)))
102 (let ((mat (la-data-to-matrix a
+mode-re
+))
103 (dp (la-data-to-vector dpars
+mode-re
+)))
106 (chol-decomp-front mat n dp
)
107 (la-matrix-to-data mat n n
+mode-re
+ result
)
108 (la-vector-to-data dp
2 +mode-re
+ dpars
))
109 (la-free-matrix mat n
)
110 (la-free-vector dp
)))
111 (list result
(second dpars
))))
116 ;;; i.e. result use by:
117 ;;; (setf (values (lu-out1 lu-out2 lu-out3)) (matlisp:lu my-matrix))
118 ;;; for solution, ...
120 ;;; (matlisp:gesv a b &opt ipivot)
124 A is a square matrix of numbers (real or complex). Computes the LU
125 decomposition of A and returns a list of the form (LU IV D FLAG), where
126 LU is a matrix with the L part in the lower triangle, the U part in the
127 upper triangle (the diagonal entries of L are taken to be 1), IV is a vector
128 describing the row permutation used, D is 1 if the number of permutations
129 is odd, -1 if even, and FLAG is T if A is numerically singular, NIL otherwise.
131 (check-square-matrix a
)
132 (let* ((n (array-dimension a
0))
133 (mode (max +mode-re
+ (la-data-mode a
)))
134 (result (list (make-array (list n n
)) (make-array n
) nil nil
)))
135 (let ((mat (la-data-to-matrix a mode
))
136 (iv (la-vector n
+mode-in
+))
137 (d (la-vector 1 +mode-re
+))
141 (setf singular
(lu-decomp-front mat n iv mode d
))
142 (la-matrix-to-data mat n n mode
(first result
))
143 (la-vector-to-data iv n
+mode-in
+ (second result
))
144 (setf (third result
) (la-get-double d
0))
145 (setf (fourth result
) (if (= singular
0.0) nil t
)))
146 (la-free-matrix mat n
)
151 (defun lu-solve (lu lb
)
153 LU is the result of (LU-DECOMP A) for a square matrix A, B is a sequence.
154 Returns the solution to the equation Ax = B. Signals an error if A is
156 (let ((la (first lu
))
158 (check-square-matrix la
)
159 (check-sequence lidx
)
162 (let* ((n (num-rows la
))
163 (result (make-sequence (if (consp lb
) 'list
'vector
) n
))
164 (a-mode (la-data-mode la
))
165 (b-mode (la-data-mode lb
)))
166 (if (/= n
(length lidx
)) (error "index sequence is wrong length"))
167 (if (/= n
(length lb
)) (error "right hand side is wrong length"))
168 (let* ((mode (max +mode-re
+ a-mode b-mode
))
169 (a (la-data-to-matrix la mode
))
170 (indx (la-data-to-vector lidx
+mode-in
+))
171 (b (la-data-to-vector lb mode
))
175 (setf singular
(lu-solve-front a n indx b mode
))
176 (la-vector-to-data b n mode result
))
178 (la-free-vector indx
)
180 (if (/= 0.0 singular
) (error "matrix is (numerically) singular"))
183 (defun determinant (a)
185 Returns the determinant of the square matrix M."
186 (let* ((lu (lu-decomp a
))
192 (flet ((fabs (x) (float (abs x
) 0.d0
)))
193 (dotimes (i n
(* d1
(exp d2
)))
195 (let* ((x (aref la i i
))
197 (if (= 0.0 magn
) (return 0.d0
))
198 (setf d1
(* d1
(/ x magn
)))
199 (setf d2
(+ d2
(log magn
))))))))
203 Returns the inverse of the the square matrix M; signals an error if M is ill
204 conditioned or singular"
205 (check-square-matrix a
)
206 (let ((n (num-rows a
))
207 (mode (max +mode-re
+ (la-data-mode a
))))
209 (let ((result (make-array (list n n
) :initial-element
0)))
212 (setf (aref result i i
) 1))
213 (let ((mat (la-data-to-matrix a mode
))
214 (inv (la-data-to-matrix result mode
))
215 (iv (la-vector n
+mode-in
+))
216 (v (la-vector n mode
))
220 (setf singular
(lu-inverse-front mat n iv v mode inv
))
221 (la-matrix-to-data inv n n mode result
))
222 (la-free-matrix mat n
)
223 (la-free-matrix inv n
)
226 (if (/= singular
0) (error "matrix is (numerically) singular"))
230 ;;;; SV Decomposition
235 A is a matrix of real numbers with at least as many rows as columns.
236 Computes the singular value decomposition of A and returns a list of the form
237 (U W V FLAG) where U and V are matrices whose columns are the left and right
238 singular vectors of A and W is the sequence of singular values of A. FLAG is T
239 if the algorithm converged, NIL otherwise."
241 (let* ((m (num-rows a
))
243 (mode (max +mode-re
+ (la-data-mode a
)))
244 (result (list (make-array (list m n
))
246 (make-array (list n n
))
248 (if (< m n
) (error "number of rows less than number of columns"))
249 (if (= mode
+mode-cx
+) (error "complex SVD not available yet"))
250 (let ((mat (la-data-to-matrix a mode
))
251 (w (la-vector n
+mode-re
+))
252 (v (la-matrix n n
+mode-re
+))
256 (setf converged
(sv-decomp-front mat m n w v
))
257 (la-matrix-to-data mat m n mode
(first result
))
258 (la-vector-to-data w n mode
(second result
))
259 (la-matrix-to-data v n n mode
(third result
))
260 (setf (fourth result
) (if (/= 0.0 converged
) t nil
)))
261 (la-free-matrix mat m
)
263 (la-free-matrix v n
))
268 ;;;; QR Decomposition
271 (defun qr-decomp (a &optional pivot
)
272 "Args: (a &optional pivot)
273 A is a matrix of real numbers with at least as many rows as columns. Computes
274 the QR factorization of A and returns the result in a list of the form (Q R).
275 If PIVOT is true the columns of X are first permuted to place insure the
276 absolute values of the diagonal elements of R are nonincreasing. In this case
277 the result includes a third element, a list of the indices of the columns in
278 the order in which they were used."
280 (let* ((m (num-rows a
))
282 (mode (max +mode-re
+ (la-data-mode a
)))
285 (list (make-array (list m n
))
286 (make-array (list n n
))
288 (list (make-array (list m n
)) (make-array (list n n
))))))
289 (if (< m n
) (error "number of rows less than number of columns"))
290 (if (= mode
+mode-cx
+) (error "complex QR decomposition not available yet"))
291 (let ((mat (la-data-to-matrix a mode
))
292 (v (la-matrix n n
+mode-re
+))
293 (jpvt (la-vector n
+mode-in
+)))
296 (qr-decomp-front mat m n v jpvt p
)
297 (la-matrix-to-data mat m n mode
(first result
))
298 (la-matrix-to-data v n n mode
(second result
))
299 (if pivot
(la-vector-to-data jpvt n
+mode-in
+ (third result
))))
300 (la-free-matrix mat m
)
302 (la-free-vector jpvt
))
306 ;;;; Estimate of Condition Number for Lower Triangular Matrix
311 Returns an estimate of the reciprocal of the L1 condition number of an upper
312 triangular matrix a."
313 (check-square-matrix a
)
314 (let ((mode (max +mode-re
+ (la-data-mode a
)))
316 (if (= mode
+mode-cx
+)
317 (error "complex condition estimate not available yet"))
318 (let ((mat (la-data-to-matrix a mode
))
321 (setf est
(rcondest-front mat n
))
322 (la-free-matrix mat n
))
326 ;;;; Make Rotation Matrix
329 (defun make-rotation (x y
&optional alpha
)
330 "Args: (x y &optional alpha)
331 Returns a rotation matrix for rotating from X to Y, or from X toward Y
332 by angle ALPHA, in radians. X and Y are sequences of the same length."
335 (if alpha
(check-one-real alpha
))
336 (let* ((n (length x
))
337 (mode (max +mode-re
+ (la-data-mode x
) (la-data-mode y
)))
338 (use-angle (if alpha
1 0))
339 (angle (if alpha
(float alpha
0.0) 0.0))
340 (result (make-array (list n n
))))
341 (if (/= n
(length y
)) (error "sequences not the same length"))
342 (if (= mode
+mode-cx
+) (error "complex data not supported yet"))
343 (let ((px (la-data-to-vector x
+mode-re
+))
344 (py (la-data-to-vector y
+mode-re
+))
345 (rot (la-matrix n n
+mode-re
+)))
348 (make-rotation-front n rot px py use-angle angle
)
349 (la-matrix-to-data rot n n
+mode-re
+ result
))
352 (la-free-matrix rot n
))
356 ;;;; Eigenvalues and Vectors
361 Returns list of list of eigenvalues and list of eigenvectors of square,
362 symmetric matrix A. Third element of result is NIL if algorithm converges.
363 If the algorithm does not converge, the third element is an integer I.
364 In this case the eigenvalues 0, ..., I are not reliable."
365 (check-square-matrix a
)
366 (let ((mode (max +mode-re
+ (la-data-mode a
)))
368 (if (= mode
+mode-cx
+) (error "matrix must be real and symmetric"))
369 (let ((evals (make-array n
))
370 (evecs (make-list (* n n
)))
371 (pa (la-data-to-vector (compound-data-seq a
) +mode-re
+))
372 (w (la-vector n
+mode-re
+))
373 (z (la-vector (* n n
) +mode-re
+))
374 (fv1 (la-vector n
+mode-re
+))
378 (setf ierr
(eigen-front pa n w z fv1
))
379 (la-vector-to-data w n
+mode-re
+ evals
)
380 (la-vector-to-data z
(* n n
) +mode-re
+ evecs
))
384 (la-free-vector fv1
))
385 (list (nreverse evals
)
386 (nreverse (mapcar #'(lambda (x) (coerce x
'vector
))
387 (split-list evecs n
)))
388 (if (/= 0 ierr
) (- n ierr
))))))
391 ;;;; Spline Interpolation
394 (defun make-smoother-args (x y xvals
)
400 (unless (integerp xvals
)
401 (check-sequence xvals
)
403 (let* ((n (length x
))
404 (ns (if (integerp xvals
) xvals
(length xvals
)))
405 (result (list (make-list ns
) (make-list ns
))))
406 (if (and y
(/= n
(length y
))) (error "sequences not the same length"))
407 (list x y n
(if (integerp xvals
) 0 1) ns xvals result
)))
409 (defun get-smoother-result (args) (seventh args
))
411 (defmacro with-smoother-data
((x y xvals is-reg
) &rest body
)
418 (unless (integerp ,xvals
)
419 (check-sequence ,xvals
)
421 (let* ((supplied (not (integerp ,xvals
)))
423 (ns (if supplied
(length ,xvals
) ,xvals
))
424 (result (list (make-list ns
) (make-list ns
))))
425 (if (and ,is-reg
(/= n
(length ,y
)))
426 (error "sequences not the same length"))
427 (if (and (not supplied
) (< ns
2))
428 (error "too few points for interpolation"))
429 (let* ((px (la-data-to-vector ,x
+mode-re
+))
430 (py (if ,is-reg
(la-data-to-vector ,y
+mode-re
+)))
432 (la-data-to-vector ,xvals
+mode-re
+)
433 (la-vector ns
+mode-re
+)))
434 (pys (la-vector ns
+mode-re
+)))
435 (unless supplied
(la-range-to-rseq n px ns pxs
))
438 (la-vector-to-data pxs ns
+mode-re
+ (first result
))
439 (la-vector-to-data pys ns
+mode-re
+ (second result
)))
441 (if ,is-reg
(la-free-vector py
))
443 (la-free-vector pys
))
446 (defun spline (x y
&key
(xvals 30))
447 "Args: (x y &key xvals)
448 Returns list of x and y values of natural cubic spline interpolation of (X,Y).
449 X must be strictly increasing. XVALS can be an integer, the number of equally
450 spaced points to use in the range of X, or it can be a sequence of points at
451 which to interpolate."
452 (with-smoother-data (x y xvals t
)
453 (let ((work (la-vector (* 2 n
) +mode-re
+))
456 (setf error
(spline-front n px py ns pxs pys work
))
457 (la-free-vector work
))
458 (if (/= error
0) (error "bad data for splines")))))
461 ;;;; Kernel Density Estimators and Smoothers
464 (defun kernel-type-code (type)
465 (cond ((eq type
'u
) 0)
470 (defun kernel-dens (x &key
(type 'b
) (width -
1.0) (xvals 30))
471 "Args: (x &key xvals width type)
472 Returns list of x and y values of kernel density estimate of X. XVALS can be an
473 integer, the number of equally spaced points to use in the range of X, or it
474 can be a sequence of points at which to interpolate. WIDTH specifies the
475 window width. TYPE specifies the lernel and should be one of the symbols G, T,
476 U or B for gaussian, triangular, uniform or bisquare. The default is B."
477 (check-one-real width
)
478 (with-smoother-data (x nil xvals nil
) ;; warning about deleting unreachable code is TRUE -- 2nd arg=nil!
479 (let ((code (kernel-type-code type
))
481 (setf error
(kernel-dens-front px n width pxs pys ns code
))
482 (if (/= 0 error
) (error "bad kernel density data")))))
484 (defun kernel-smooth (x y
&key
(type 'b
) (width -
1.0) (xvals 30))
485 "Args: (x y &key xvals width type)
486 Returns list of x and y values of kernel smooth of (X,Y). XVALS can be an
487 integer, the number of equally spaced points to use in the range of X, or it
488 can be a sequence of points at which to interpolate. WIDTH specifies the
489 window width. TYPE specifies the lernel and should be one of the symbols G, T,
490 U or B for Gaussian, triangular, uniform or bisquare. The default is B."
491 (check-one-real width
)
492 (with-smoother-data (x y xvals t
)
493 (let ((code (kernel-type-code type
))
495 (kernel-smooth-front px py n width pxs pys ns code
)
496 ;; if we get the Lisp version ported from C, uncomment below and
497 ;; comment above. (thanks to Carlos Ungil for the initial CFFI
499 ;;(kernel-smooth-Cport px py n width pxs pys ns code)
500 (if (/= 0 error
) (error "bad kernel density data")))))
504 (defun kernel-smooth-Cport (px py n width
;;wts wds ;; see above for mismatch?
506 "Port of kernel_smooth (Lib/kernel.c) to Lisp.
507 FIXME:kernel-smooth-Cport : This is broken.
508 Until this is fixed, we are using Luke's C code and CFFI as glue."
510 ((and (< n
2) (<= width
0)) 1.0)
511 (t (let* ((xmin (min px
))
513 (width (/ (- xmax xmin
) (+ 1.0 (log n
)))))
514 (dotimes (i (- ns
1))
518 (dotimes (j (- n
1)) )
519 ;;;possible nasty errors...
522 ((lwidth (if wds
(* width
(aref wds j
)) width
))
523 (lwt (* (kernel-Cport (aref xs i
) (aref px j
) lwidth ktype
) ;; px?
524 (if wts
(aref wts j
) 1.0))))
525 (setf wsum
(+ wsum lwt
))
526 (setf ysum
(if py
(+ ysum
(* lwt
(aref py j
)))))) ;; py? y?
536 (defun kernel-Cport (x y w ktype
)
537 "Port of kernel() (Lib/kernel.c) to Lisp.
538 x,y,w are doubles, type is an integer"
542 (cond ((eq ktype
"B")
547 (/ (/ (* 15.0 (* (- 1.0 (* 4 z z
)) ;; k/w
548 (- 1.0 (* 4 z z
)))) ;; k/w
553 (let* ((w (* w
0.25))
555 (k (/ (exp (* -
0.5 z z
))
561 (k (if (< (abs z
) 0.5)
566 (cond ((and (> z -
1.0)
577 ;;;; Lowess Smoother Interface
580 (defun |base-lowess|
(s1 s2 f nsteps delta
)
586 (check-one-fixnum nsteps
)
587 (check-one-real delta
)
588 (let* ((n (length s1
))
589 (result (make-list n
)))
590 (if (/= n
(length s2
)) (error "sequences not the same length"))
591 (let ((x (la-data-to-vector s1
+mode-re
+))
592 (y (la-data-to-vector s2
+mode-re
+))
593 (ys (la-vector n
+mode-re
+))
594 (rw (la-vector n
+mode-re
+))
595 (res (la-vector n
+mode-re
+))
599 (setf error
(base-lowess-front x y n f nsteps delta ys rw res
))
600 (la-vector-to-data ys n
+mode-re
+ result
))
605 (la-free-vector res
))
606 (if (/= error
0) (error "bad data for lowess"))
610 static LVAL add_contour_point
(i, j
, k
, l
, x
, y
, z
, v
, result
)
620 if
((z[i][j] <= v && v < z[k][l]) || (z[k][l] <= v && v < z[i][j])) {
623 p = (v - z[i][j]) / (z[k][l] - z[i][j]);
625 rplaca(pt, cvflonum((FLOTYPE) (q * x[i] + p * x[k])));
626 rplaca
(cdr(pt), cvflonum
((FLOTYPE) (q * y
[j] + p * y[l])));
627 result = cons(pt, result);
633 LVAL xssurface_contour()
635 LVAL s1, s2, mat, result;
641 s1 = xsgetsequence();
642 s2 = xsgetsequence();
644 v = makedouble(xlgetarg());
647 n = seqlen(s1); m = seqlen(s2);
648 if (n != numrows(mat) || m != numcols(mat)) xlfail("dimensions do not match");
649 if (data_mode(s1) == CX || data_mode(s2) == CX || data_mode(mat) == CX)
650 xlfail("data must be real");
652 x = (RVector) data_to_vector(s1, RE);
653 y = (RVector) data_to_vector(s2, RE);
654 z = (RMatrix) data_to_matrix(mat, RE);
658 for (i = 0; i < n - 1; i++) {
659 for (j = 0; j < m - 1; j++) {
660 result = add_contour_point(i, j, i, j+1, x, y, z, v, result);
661 result = add_contour_point(i, j+1, i+1, j+1, x, y, z, v, result);
662 result = add_contour_point(i+1, j+1, i+1, j, x, y, z, v, result);
663 result = add_contour_point(i+1, j, i, j, x, y, z, v, result);
680 ;;; ??replace with matlisp:fft and matlisp:ifft (the latter for inverse mapping)
682 (defun fft (x &optional inverse)
683 "Args: (x &optional inverse)
684 Returns unnormalized Fourier transform of X, or inverse transform if INVERSE
687 (let* ((n (length x))
688 (mode (la-data-mode x))
689 (isign (if inverse -1 1))
690 (result (if (consp x) (make-list n) (make-array n))))
691 (let ((px (la-data-to-vector x +mode-cx+))
692 (work (la-vector (+ (* 4 n) 15) +mode-re+)))
695 (fft-front n px work isign)
696 (la-vector-to-data px n +mode-cx+ result))
698 (la-free-vector work))
702 ;;; SWEEP Operator: FIXME: use matlisp
705 (defun make-sweep-front (x y w n p mode has_w x_mean result)
706 (declare (fixnum n p mode has_w))
721 (has-w (if (/= 0 has_w) t nil))
723 (declare (long-float val dxi dyi dv dw sum_w dxik dxjk dyj
724 dx_meani dx_meanj dy_mean))
725 ;; (declare-double val dxi dyi dv dw sum_w dxik dxjk dyj
726 ;; dx_meani dx_meanj dy_mean)
728 (if (> mode RE) (error "not supported for complex data yet"))
730 (setf x_data (compound-data-seq x))
731 (setf result_data (compound-data-seq result))
733 ;; find the mean of y
738 (setf dyi (makedouble (aref y i)))
740 (setf dw (makedouble (aref w i)))
742 (setf dyi (* dyi dw)))
744 (if (not has-w) (setf sum_w (float n 0.0)))
745 (if (<= sum_w 0.0) (error "non positive sum of weights"))
746 (setf dy_mean (/ val sum_w))
748 ;; find the column means
754 (setf dxi (makedouble (aref x_data (+ (* p i) j))))
756 (setf dw (makedouble (aref w i)))
757 (setf dxi (* dxi dw)))
759 (setf (aref x_mean j) (/ val sum_w)))
761 ;; put 1/sum_w in topleft, means on left, minus means on top
762 (setf (aref result_data 0) (/ 1.0 sum_w))
765 (setf dxi (makedouble (aref x_mean i)))
766 (setf (aref result_data (+ i 1)) (- dxi))
767 (setf (aref result_data (* (+ i 1) (+ p 2))) dxi))
768 (setf (aref result_data (+ p 1)) (- dy_mean))
769 (setf (aref result_data (* (+ p 1) (+ p 2))) dy_mean)
771 ;; put sums of adjusted cross products in body
779 (setf dxik (makedouble (aref x_data (+ (* p k) i))))
780 (setf dxjk (makedouble (aref x_data (+ (* p k) j))))
781 (setf dx_meani (makedouble (aref x_mean i)))
782 (setf dx_meanj (makedouble (aref x_mean j)))
783 (setf dv (* (- dxik dx_meani) (- dxjk dx_meanj)))
785 (setf dw (makedouble (aref w k)))
788 (setf (aref result_data (+ (* (+ i 1) (+ p 2)) (+ j 1))) val)
789 (setf (aref result_data (+ (* (+ j 1) (+ p 2)) (+ i 1))) val))
793 (setf dxik (makedouble (aref x_data (+ (* p j) i))))
794 (setf dyj (makedouble (aref y j)))
795 (setf dx_meani (makedouble (aref x_mean i)))
796 (setf dv (* (- dxik dx_meani) (- dyj dy_mean)))
798 (setf dw (makedouble (aref w j)))
801 (setf (aref result_data (+ (* (+ i 1) (+ p 2)) (+ p 1))) val)
802 (setf (aref result_data (+ (* (+ p 1) (+ p 2)) (+ i 1))) val))
806 (setf dyj (makedouble (aref y j)))
807 (setf dv (* (- dyj dy_mean) (- dyj dy_mean)))
809 (setf dw (makedouble (aref w j)))
812 (setf (aref result_data (+ (* (+ p 1) (+ p 2)) (+ p 1))) val)))
814 ;;; FIXME: use matlisp
815 (defun sweep-in-place-front (a rows cols mode k tol)
816 "Sweep algorithm for linear regression."
817 (declare (long-float tol))
818 (declare (fixnum rows cols mode k))
826 (declare (long-float pivot aij aik akj akk))
828 (if (> mode RE) (error "not supported for complex data yet"))
829 (if (or (< k 0) (>= k rows) (>= k cols)) (error "index out of range"))
831 (setf tol (max tol machine-epsilon))
832 (setf data (compound-data-seq a))
834 (setf pivot (makedouble (aref data (+ (* cols k) k))))
837 ((or (> pivot tol) (< pivot (- tol)))
842 (when (and (/= i k) (/= j k))
843 (setf aij (makedouble (aref data (+ (* cols i) j))))
844 (setf aik (makedouble (aref data (+ (* cols i) k))))
845 (setf akj (makedouble (aref data (+ (* cols k) j))))
846 (setf aij (- aij (/ (* aik akj) pivot)))
847 (setf (aref data (+ (* cols i) j)) aij))))
851 (setf aik (makedouble (aref data (+ (* cols i) k))))
853 (setf aik (/ aik pivot))
854 (setf (aref data (+ (* cols i) k)) aik)))
858 (setf akj (makedouble (aref data (+ (* cols k) j))))
860 (setf akj (- (/ akj pivot)))
861 (setf (aref data (+ (* cols k) j)) akj)))
863 (setf akk (/ 1.0 pivot))
864 (setf (aref data (+ (* cols k) k)) akk)
868 ;; FIXME: use matlisp
869 (defun make-sweep-matrix (x y &optional w)
870 "Args: (x y &optional weights)
871 X is matrix, Y and WEIGHTS are sequences. Returns the sweep matrix of the
872 (weighted) regression of Y on X"
875 (if w (check-sequence w))
876 (let ((n (num-rows x))
878 (if (/= n (length y)) (error "dimensions do not match"))
879 (if (and w (/= n (length w))) (error "dimensions do not match"))
880 (let ((mode (max (la-data-mode x)
882 (if w (la-data-mode w) 0)))
883 (result (make-array (list (+ p 2) (+ p 2))))
884 (x-mean (make-array p))
885 (y (coerce y 'vector))
886 (w (if w (coerce w 'vector)))
888 (make-sweep-front x y w n p mode has-w x-mean result)
891 (defun sweep-in-place (a k tol)
895 (let ((rows (num-rows a))
897 (mode (la-data-mode a)))
898 (let ((swept (sweep-in-place-front a rows cols mode k tol)))
899 (if (/= 0 swept) t nil))))
901 (defun sweep-operator (a columns &optional tolerances)
902 "Args: (a indices &optional tolerances)
903 A is a matrix, INDICES a sequence of the column indices to be swept. Returns
904 a list of the swept result and the list of the columns actually swept. (See
905 MULTREG documentation.) If supplied, TOLERANCES should be a list of real
906 numbers the same length as INDICES. An index will only be swept if its pivot
907 element is larger than the corresponding element of TOLERANCES."
909 (check-sequence columns)
910 (if tolerances (check-sequence tolerances))
912 (check-fixnum columns)
913 (if tolerances (check-real tolerances))
915 (result (copy-array a))
917 (columns (coerce columns 'list) (cdr columns))
918 (tolerances (if (consp tolerances) (coerce tolerances 'list))
919 (if (consp tolerances) (cdr tolerances))))
920 ((null columns) (list result swept-columns))
921 (let ((col (first columns))
922 (tol (if (consp tolerances) (first tolerances) tol)))
923 (if (sweep-in-place result col tol)
924 (setf swept-columns (cons col swept-columns))))))
933 (defun ax+y (a x y &optional lower)
934 "Args (a x y &optional lower)
935 Returns (+ (matmult A X) Y). If LOWER is not nil, A is taken to be lower
937 This could probably be made more efficient."
938 (check-square-matrix a)
944 (let* ((n (num-rows a))
945 (result (make-list n))
946 (a (compound-data-seq a)))
948 (if (or (/= n (length x)) (/= n (length y)))
949 (error "dimensions do not match"))
950 (do* ((tx (make-next-element x) (make-next-element x))
951 (ty (make-next-element y))
952 (tr (make-next-element result))
954 (start 0 (+ start n))
955 (end (if lower (+ i 1) n) (if lower (+ i 1) n)))
957 (declare (fixnum i start end))
958 (let ((val (get-next-element ty i)))
961 (setf val (+ val (* (get-next-element tx j)
962 (aref a (+ start j))))))
963 (set-next-element tr i val)))))