clean up targets and dependencies; cleaner CL despite sideeffect functions
[CommonLispStat.git] / linalg.lsp
blobbce7e804b14d078571a301d96bf2ef3322f8ffae
1 ;;; -*- mode: lisp -*-
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.
5 ;;;
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.
14 ;;;
15 ;;; #3 - Use a lisp-based matrix system drop-in? (matlisp, femlisp, clem, ...?)
16 ;;;
19 ;;;; linalg -- Lisp-Stat interface to basic linear algebra routines.
20 ;;;;
21 ;;;; Copyright (c) 1991, by Luke Tierney. Permission is granted for
22 ;;;; unrestricted use.
24 ;;;;
25 ;;;; Package Setup
26 ;;;;
28 (defpackage :lisp-stat-linalg
29 (:use :common-lisp
30 :lisp-stat-math
31 :lisp-stat-types
32 :lisp-stat-matrix)
33 (:shadowing-import-from :lisp-stat-math
34 expt + - * / ** mod rem abs 1+ 1- log exp sqrt sin cos tan
35 asin acos atan sinh cosh tanh asinh acosh atanh float random
36 truncate floor ceiling round minusp zerop plusp evenp oddp
37 < <= = /= >= > complex conjugate realpart imagpart phase
38 min max logand logior logxor lognot ffloor fceiling
39 ftruncate fround signum cis)
40 (:export chol-decomp lu-decomp lu-solve determinant inverse sv-decomp
41 qr-decomp rcondest make-rotation spline kernel-dens kernel-smooth
42 fft make-sweep-matrix sweep-operator ax+y eigen))
44 (in-package #:lisp-stat-linalg)
46 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
47 ;;;;
48 ;;;; Lisp to C number conversion and checking
49 ;;;;
50 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
52 ;;;;
53 ;;;; Lisp to/from C sequence and matrix conversion and checking
54 ;;;;
56 (defun is-cons (a)
57 "FIXME:AJR this is not used anywhere?"
58 (if (consp a) 1 0))
60 (defun check-fixnum (a)
61 (if (/= 0 (la-data-mode a)) (error "not an integer sequence - ~s" a)))
63 (defun check-real (data)
64 (let ((data (compound-data-seq data)))
65 (cond
66 ((vectorp data)
67 (let ((n (length data)))
68 (declare (fixnum n))
69 (dotimes (i n)
70 (declare (fixnum i))
71 (check-one-real (aref data i)))))
72 ((consp data) (dolist (x data) (check-one-real x)))
73 (t (error "bad sequence - ~s" data)))))
75 (defun vec-assign (a i x) (setf (aref a i) x))
77 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
78 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
79 ;;;;
80 ;;;; Lisp Interfaces to Linear Algebra Routines
81 ;;;;
82 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
83 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
85 ;;; FIXME: use dpbt[f2|rf], dpbstf, dpot[f2|rf]; dpptrf, zpbstf, zpbt[f2|rf]
86 ;;; remember: factorization = decomposition, depending on training.
88 (defun chol-decomp (a &optional (maxoffl 0.0))
89 "Args: (a)
90 Modified Cholesky decomposition. A should be a square, symmetric matrix.
91 Computes lower triangular matrix L such that L L^T = A + D where D is a
92 diagonal matrix. If A is strictly positive definite D will be zero.
93 Otherwise, D is as small as possible to make A + D numerically strictly
94 positive definite. Returns a list (L (max D))."
95 (check-square-matrix a)
96 (check-real a)
97 (let* ((n (array-dimension a 0))
98 (result (make-array (list n n)))
99 (dpars (list maxoffl 0.0)))
100 (check-real dpars)
101 (let ((mat (la-data-to-matrix a mode-re))
102 (dp (la-data-to-vector dpars mode-re)))
103 (unwind-protect
104 (progn
105 (chol-decomp-front mat n dp)
106 (la-matrix-to-data mat n n mode-re result)
107 (la-vector-to-data dp 2 mode-re dpars))
108 (la-free-matrix mat n)
109 (la-free-vector dp)))
110 (list result (second dpars))))
113 ;;; REPLACE with
114 ;;; (matlisp:lu M)
115 ;;; i.e. result use by:
116 ;;; (setf (values (lu-out1 lu-out2 lu-out3)) (matlisp:lu my-matrix))
117 ;;; for solution, ...
118 ;;; for lu-solve:
119 ;;; (matlisp:gesv a b &opt ipivot)
121 (defun lu-decomp (a)
122 "Args: (a)
123 A is a square matrix of numbers (real or complex). Computes the LU
124 decomposition of A and returns a list of the form (LU IV D FLAG), where
125 LU is a matrix with the L part in the lower triangle, the U part in the
126 upper triangle (the diagonal entries of L are taken to be 1), IV is a vector
127 describing the row permutation used, D is 1 if the number of permutations
128 is odd, -1 if even, and FLAG is T if A is numerically singular, NIL otherwise.
129 Used bu LU-SOLVE."
130 (check-square-matrix a)
131 (let* ((n (array-dimension a 0))
132 (mode (max mode-re (la-data-mode a)))
133 (result (list (make-array (list n n)) (make-array n) nil nil)))
134 (let ((mat (la-data-to-matrix a mode))
135 (iv (la-vector n mode-in))
136 (d (la-vector 1 mode-re))
137 (singular 0))
138 (unwind-protect
139 (progn
140 (setf singular (lu-decomp-front mat n iv mode d))
141 (la-matrix-to-data mat n n mode (first result))
142 (la-vector-to-data iv n mode-in (second result))
143 (setf (third result) (la-get-double d 0))
144 (setf (fourth result) (if (= singular 0.0) nil t)))
145 (la-free-matrix mat n)
146 (la-free-vector iv)
147 (la-free-vector d)))
148 result))
150 (defun lu-solve (lu lb)
151 "Args: (lu b)
152 LU is the result of (LU-DECOMP A) for a square matrix A, B is a sequence.
153 Returns the solution to the equation Ax = B. Signals an error if A is
154 singular."
155 (let ((la (first lu))
156 (lidx (second lu)))
157 (check-square-matrix la)
158 (check-sequence lidx)
159 (check-sequence lb)
160 (check-fixnum lidx)
161 (let* ((n (num-rows la))
162 (result (make-sequence (if (consp lb) 'list 'vector) n))
163 (a-mode (la-data-mode la))
164 (b-mode (la-data-mode lb)))
165 (if (/= n (length lidx)) (error "index sequence is wrong length"))
166 (if (/= n (length lb)) (error "right hand side is wrong length"))
167 (let* ((mode (max mode-re a-mode b-mode))
168 (a (la-data-to-matrix la mode))
169 (indx (la-data-to-vector lidx mode-in))
170 (b (la-data-to-vector lb mode))
171 (singular 0))
172 (unwind-protect
173 (progn
174 (setf singular (lu-solve-front a n indx b mode))
175 (la-vector-to-data b n mode result))
176 (la-free-matrix a n)
177 (la-free-vector indx)
178 (la-free-vector b))
179 (if (/= 0.0 singular) (error "matrix is (numerically) singular"))
180 result))))
182 (defun determinant (a)
183 "Args: (m)
184 Returns the determinant of the square matrix M."
185 (let* ((lu (lu-decomp a))
186 (la (first lu))
187 (n (num-rows a))
188 (d1 (third lu))
189 (d2 0.d0))
190 (declare (fixnum n))
191 (flet ((fabs (x) (float (abs x) 0.d0)))
192 (dotimes (i n (* d1 (exp d2)))
193 (declare (fixnum i))
194 (let* ((x (aref la i i))
195 (magn (fabs x)))
196 (if (= 0.0 magn) (return 0.d0))
197 (setf d1 (* d1 (/ x magn)))
198 (setf d2 (+ d2 (log magn))))))))
200 (defun inverse (a)
201 "Args: (m)
202 Returns the inverse of the the square matrix M; signals an error if M is ill
203 conditioned or singular"
204 (check-square-matrix a)
205 (let ((n (num-rows a))
206 (mode (max mode-re (la-data-mode a))))
207 (declare (fixnum n))
208 (let ((result (make-array (list n n) :initial-element 0)))
209 (dotimes (i n)
210 (declare (fixnum i))
211 (setf (aref result i i) 1))
212 (let ((mat (la-data-to-matrix a mode))
213 (inv (la-data-to-matrix result mode))
214 (iv (la-vector n mode-in))
215 (v (la-vector n mode))
216 (singular 0))
217 (unwind-protect
218 (progn
219 (setf singular (lu-inverse-front mat n iv v mode inv))
220 (la-matrix-to-data inv n n mode result))
221 (la-free-matrix mat n)
222 (la-free-matrix inv n)
223 (la-free-vector iv)
224 (la-free-vector v))
225 (if (/= singular 0) (error "matrix is (numerically) singular"))
226 result))))
228 ;;;;
229 ;;;; SV Decomposition
230 ;;;;
232 (defun sv-decomp (a)
233 "Args: (a)
234 A is a matrix of real numbers with at least as many rows as columns.
235 Computes the singular value decomposition of A and returns a list of the form
236 (U W V FLAG) where U and V are matrices whose columns are the left and right
237 singular vectors of A and W is the sequence of singular values of A. FLAG is T
238 if the algorithm converged, NIL otherwise."
239 (check-matrix a)
240 (let* ((m (num-rows a))
241 (n (num-cols a))
242 (mode (max mode-re (la-data-mode a)))
243 (result (list (make-array (list m n))
244 (make-array n)
245 (make-array (list n n))
246 nil)))
247 (if (< m n) (error "number of rows less than number of columns"))
248 (if (= mode mode-cx) (error "complex SVD not available yet"))
249 (let ((mat (la-data-to-matrix a mode))
250 (w (la-vector n mode-re))
251 (v (la-matrix n n mode-re))
252 (converged 0))
253 (unwind-protect
254 (progn
255 (setf converged (sv-decomp-front mat m n w v))
256 (la-matrix-to-data mat m n mode (first result))
257 (la-vector-to-data w n mode (second result))
258 (la-matrix-to-data v n n mode (third result))
259 (setf (fourth result) (if (/= 0.0 converged) t nil)))
260 (la-free-matrix mat m)
261 (la-free-vector w)
262 (la-free-matrix v n))
263 result)))
266 ;;;;
267 ;;;; QR Decomposition
268 ;;;;
270 (defun qr-decomp (a &optional pivot)
271 "Args: (a &optional pivot)
272 A is a matrix of real numbers with at least as many rows as columns. Computes
273 the QR factorization of A and returns the result in a list of the form (Q R).
274 If PIVOT is true the columns of X are first permuted to place insure the
275 absolute values of the diagonal elements of R are nonincreasing. In this case
276 the result includes a third element, a list of the indices of the columns in
277 the order in which they were used."
278 (check-matrix a)
279 (let* ((m (num-rows a))
280 (n (num-cols a))
281 (mode (max mode-re (la-data-mode a)))
282 (p (if pivot 1 0))
283 (result (if pivot
284 (list (make-array (list m n))
285 (make-array (list n n))
286 (make-array n))
287 (list (make-array (list m n)) (make-array (list n n))))))
288 (if (< m n) (error "number of rows less than number of columns"))
289 (if (= mode mode-cx) (error "complex QR decomposition not available yet"))
290 (let ((mat (la-data-to-matrix a mode))
291 (v (la-matrix n n mode-re))
292 (jpvt (la-vector n mode-in)))
293 (unwind-protect
294 (progn
295 (qr-decomp-front mat m n v jpvt p)
296 (la-matrix-to-data mat m n mode (first result))
297 (la-matrix-to-data v n n mode (second result))
298 (if pivot (la-vector-to-data jpvt n mode-in (third result))))
299 (la-free-matrix mat m)
300 (la-free-matrix v n)
301 (la-free-vector jpvt))
302 result)))
304 ;;;;
305 ;;;; Estimate of Condition Number for Lower Triangular Matrix
306 ;;;;
308 (defun rcondest (a)
309 "Args: (a)
310 Returns an estimate of the reciprocal of the L1 condition number of an upper
311 triangular matrix a."
312 (check-square-matrix a)
313 (let ((mode (max mode-re (la-data-mode a)))
314 (n (num-rows a)))
315 (if (= mode mode-cx)
316 (error "complex condition estimate not available yet"))
317 (let ((mat (la-data-to-matrix a mode))
318 (est 0.0))
319 (unwind-protect
320 (setf est (rcondest-front mat n))
321 (la-free-matrix mat n))
322 est)))
324 ;;;;
325 ;;;; Make Rotation Matrix
326 ;;;;
328 (defun make-rotation (x y &optional alpha)
329 "Args: (x y &optional alpha)
330 Returns a rotation matrix for rotating from X to Y, or from X toward Y
331 by angle ALPHA, in radians. X and Y are sequences of the same length."
332 (check-sequence x)
333 (check-sequence y)
334 (if alpha (check-one-real alpha))
335 (let* ((n (length x))
336 (mode (max mode-re (la-data-mode x) (la-data-mode y)))
337 (use-angle (if alpha 1 0))
338 (angle (if alpha (float alpha 0.0) 0.0))
339 (result (make-array (list n n))))
340 (if (/= n (length y)) (error "sequences not the same length"))
341 (if (= mode mode-cx) (error "complex data not supported yet"))
342 (let ((px (la-data-to-vector x mode-re))
343 (py (la-data-to-vector y mode-re))
344 (rot (la-matrix n n mode-re)))
345 (unwind-protect
346 (progn
347 (make-rotation-front n rot px py use-angle angle)
348 (la-matrix-to-data rot n n mode-re result))
349 (la-free-vector px)
350 (la-free-vector py)
351 (la-free-matrix rot n))
352 result)))
354 ;;;;
355 ;;;; Eigenvalues and Vectors
356 ;;;;
358 (defun eigen (a)
359 "Args: (a)
360 Returns list of list of eigenvalues and list of eigenvectors of square,
361 symmetric matrix A. Third element of result is NIL if algorithm converges.
362 If the algorithm does not converge, the third element is an integer I.
363 In this case the eigenvalues 0, ..., I are not reliable."
364 (check-square-matrix a)
365 (let ((mode (max mode-re (la-data-mode a)))
366 (n (num-rows a)))
367 (if (= mode mode-cx) (error "matrix must be real and symmetric"))
368 (let ((evals (make-array n))
369 (evecs (make-list (* n n)))
370 (pa (la-data-to-vector (compound-data-seq a) mode-re))
371 (w (la-vector n mode-re))
372 (z (la-vector (* n n) mode-re))
373 (fv1 (la-vector n mode-re))
374 (ierr 0))
375 (unwind-protect
376 (progn
377 (setf ierr (eigen-front pa n w z fv1))
378 (la-vector-to-data w n mode-re evals)
379 (la-vector-to-data z (* n n) mode-re evecs))
380 (la-free-vector pa)
381 (la-free-vector z)
382 (la-free-vector w)
383 (la-free-vector fv1))
384 (list (nreverse evals)
385 (nreverse (mapcar #'(lambda (x) (coerce x 'vector))
386 (split-list evecs n)))
387 (if (/= 0 ierr) (- n ierr))))))
389 ;;;;
390 ;;;; Spline Interpolation
391 ;;;;
393 (defun make-smoother-args (x y xvals)
394 (check-sequence x)
395 (check-real x)
396 (when y
397 (check-sequence y)
398 (check-real y))
399 (unless (integerp xvals)
400 (check-sequence xvals)
401 (check-real xvals))
402 (let* ((n (length x))
403 (ns (if (integerp xvals) xvals (length xvals)))
404 (result (list (make-list ns) (make-list ns))))
405 (if (and y (/= n (length y))) (error "sequences not the same length"))
406 (list x y n (if (integerp xvals) 0 1) ns xvals result)))
408 (defun get-smoother-result (args) (seventh args))
410 (defmacro with-smoother-data ((x y xvals is-reg) &rest body)
411 `(progn
412 (check-sequence ,x)
413 (check-real ,x)
414 (when ,is-reg
415 (check-sequence ,y)
416 (check-real ,y))
417 (unless (integerp ,xvals)
418 (check-sequence ,xvals)
419 (check-real ,xvals))
420 (let* ((supplied (not (integerp ,xvals)))
421 (n (length ,x))
422 (ns (if supplied (length ,xvals) ,xvals))
423 (result (list (make-list ns) (make-list ns))))
424 (if (and ,is-reg (/= n (length ,y)))
425 (error "sequences not the same length"))
426 (if (and (not supplied) (< ns 2))
427 (error "too few points for interpolation"))
428 (let* ((px (la-data-to-vector ,x mode-re))
429 (py (if ,is-reg (la-data-to-vector ,y mode-re)))
430 (pxs (if supplied
431 (la-data-to-vector ,xvals mode-re)
432 (la-vector ns mode-re)))
433 (pys (la-vector ns mode-re)))
434 (unless supplied (la-range-to-rseq n px ns pxs))
435 (unwind-protect
436 (progn ,@body
437 (la-vector-to-data pxs ns mode-re (first result))
438 (la-vector-to-data pys ns mode-re (second result)))
439 (la-free-vector px)
440 (if ,is-reg (la-free-vector py))
441 (la-free-vector pxs)
442 (la-free-vector pys))
443 result))))
445 (defun spline (x y &key (xvals 30))
446 "Args: (x y &key xvals)
447 Returns list of x and y values of natural cubic spline interpolation of (X,Y).
448 X must be strictly increasing. XVALS can be an integer, the number of equally
449 spaced points to use in the range of X, or it can be a sequence of points at
450 which to interpolate."
451 (with-smoother-data (x y xvals t)
452 (let ((work (la-vector (* 2 n) mode-re))
453 (error 0))
454 (unwind-protect
455 (setf error (spline-front n px py ns pxs pys work))
456 (la-free-vector work))
457 (if (/= error 0) (error "bad data for splines")))))
459 ;;;;
460 ;;;; Kernel Density Estimators and Smoothers
461 ;;;;
463 (defun kernel-type-code (type)
464 (cond ((eq type 'u) 0)
465 ((eq type 't) 1)
466 ((eq type 'g) 2)
467 (t 3)))
469 (defun kernel-dens (x &key (type 'b) (width -1.0) (xvals 30))
470 "Args: (x &key xvals width type)
471 Returns list of x and y values of kernel density estimate of X. XVALS can be an
472 integer, the number of equally spaced points to use in the range of X, or it
473 can be a sequence of points at which to interpolate. WIDTH specifies the
474 window width. TYPE specifies the lernel and should be one of the symbols G, T,
475 U or B for gaussian, triangular, uniform or bisquare. The default is B."
476 (check-one-real width)
477 (with-smoother-data (x nil xvals nil) ;; warning about deleting unreachable code is TRUE -- 2nd arg=nil!
478 (let ((code (kernel-type-code type))
479 (error 0))
480 (setf error (kernel-dens-front px n width pxs pys ns code))
481 (if (/= 0 error) (error "bad kernel density data")))))
483 (defun kernel-smooth (x y &key (type 'b) (width -1.0) (xvals 30))
484 "Args: (x y &key xvals width type)
485 Returns list of x and y values of kernel smooth of (X,Y). XVALS can be an
486 integer, the number of equally spaced points to use in the range of X, or it
487 can be a sequence of points at which to interpolate. WIDTH specifies the
488 window width. TYPE specifies the lernel and should be one of the symbols G, T,
489 U or B for Gaussian, triangular, uniform or bisquare. The default is B."
490 (check-one-real width)
491 (with-smoother-data (x y xvals t)
492 (let ((code (kernel-type-code type))
493 (error 0))
494 ;;(kernel-smooth-front px py n width pxs pys ns code)
495 (kernel-smooth-Cport px py n width pxs pys ns code)
496 (if (/= 0 error) (error "bad kernel density data")))))
500 (defun kernel-smooth-Cport (px py n width ;;wts wds ;; see above for mismatch?
501 xs ys ns ktype)
502 "Port of kernel_smooth (Lib/kernel.c) to Lisp.
503 FIXME:kernel-smooth-Cport"
504 (cond ((< n 1) 1.0)
505 ((and (< n 2) (<= width 0)) 1.0)
506 (t (let* ((xmin (min px))
507 (xmax (max px))
508 (width (/ (- xmax xmin) (+ 1.0 (log n)))))
509 (dotimes (i (- ns 1))
510 (setf (aref ys i)
511 (let ((wsum 0.0)
512 (ysum 0.0))
513 (dotimes (j (- n 1)) )
514 ;;;possible nasty errors...
516 (let*
517 ((lwidth (if wds (* width (aref wds j)) width))
518 (lwt (* (kernel-Cport (aref xs i) (aref px j) lwidth ktype) ;; px?
519 (if wts (aref wts j) 1.0))))
520 (setf wsum (+ wsum lwt))
521 (setf ysum (if py (+ ysum (* lwt (aref py j)))))) ;; py? y?
523 ;;; end of errors
524 (if py
525 (if (> wsum 0.0)
526 (/ ysum wsum)
527 0.0)
528 (/ wsum n)))))
529 (values ys)))))
531 (defun kernel-Cport (x y w ktype)
532 "Port of kernel() (Lib/kernel.c) to Lisp.
533 x,y,w are doubles, type is an integer"
534 (if (<= w 0.0)
536 (let ((z (- x y)))
537 (cond ((eq ktype "B")
538 (let* ((w (* w 2.0))
539 (z (* z 0.5)))
540 (if (and (> z -0.5)
541 (< z 0.5))
542 (/ (/ (* 15.0 (* (- 1.0 (* 4 z z)) ;; k/w
543 (- 1.0 (* 4 z z)))) ;; k/w
544 8.0)
546 0)))
547 ((eq ktype "G")
548 (let* ((w (* w 0.25))
549 (z (* z 4.0))
550 (k (/ (exp (* -0.5 z z))
551 (sqrt (* 2 PI)))))
552 (/ k w)))
553 ((eq ktype "U")
554 (let* ((w (* 1.5 w))
555 (z (* z 0.75))
556 (k (if (< (abs z) 0.5)
558 0.0)))
559 (/ k w)))
560 ((eq ktype "T")
561 (cond ((and (> z -1.0)
562 (< z 0.0))
563 (+ 1.0 z)) ;; k
564 ((and (> z 0.0)
565 (< z 1.0))
566 (- 1.0 z)) ;; k
567 (t 0.0)))
568 (t (values 0.0))))))
571 ;;;;
572 ;;;; Lowess Smoother Interface
573 ;;;;
575 (defun |base-lowess| (s1 s2 f nsteps delta)
576 (check-sequence s1)
577 (check-sequence s2)
578 (check-real s1)
579 (check-real s2)
580 (check-one-real f)
581 (check-one-fixnum nsteps)
582 (check-one-real delta)
583 (let* ((n (length s1))
584 (result (make-list n)))
585 (if (/= n (length s2)) (error "sequences not the same length"))
586 (let ((x (la-data-to-vector s1 mode-re))
587 (y (la-data-to-vector s2 mode-re))
588 (ys (la-vector n mode-re))
589 (rw (la-vector n mode-re))
590 (res (la-vector n mode-re))
591 (error 0))
592 (unwind-protect
593 (progn
594 (setf error (base-lowess-front x y n f nsteps delta ys rw res))
595 (la-vector-to-data ys n mode-re result))
596 (la-free-vector x)
597 (la-free-vector y)
598 (la-free-vector ys)
599 (la-free-vector rw)
600 (la-free-vector res))
601 (if (/= error 0) (error "bad data for lowess"))
602 result)))
605 static LVAL add_contour_point(i, j, k, l, x, y, z, v, result)
606 int i, j, k, l;
607 RVector x, y;
608 RMatrix z;
609 double v;
610 LVAL result;
612 LVAL pt;
613 double p, q;
615 if ((z[i][j] <= v && v < z[k][l]) || (z[k][l] <= v && v < z[i][j])) {
616 xlsave(pt);
617 pt = mklist(2, NIL);
618 p = (v - z[i][j]) / (z[k][l] - z[i][j]);
619 q = 1.0 - p;
620 rplaca(pt, cvflonum((FLOTYPE) (q * x[i] + p * x[k])));
621 rplaca(cdr(pt), cvflonum((FLOTYPE) (q * y[j] + p * y[l])));
622 result = cons(pt, result);
623 xlpop();
625 return(result);
628 LVAL xssurface_contour()
630 LVAL s1, s2, mat, result;
631 RVector x, y;
632 RMatrix z;
633 double v;
634 int i, j, n, m;
636 s1 = xsgetsequence();
637 s2 = xsgetsequence();
638 mat = xsgetmatrix();
639 v = makedouble(xlgetarg());
640 xllastarg();
642 n = seqlen(s1); m = seqlen(s2);
643 if (n != numrows(mat) || m != numcols(mat)) xlfail("dimensions do not match");
644 if (data_mode(s1) == CX || data_mode(s2) == CX || data_mode(mat) == CX)
645 xlfail("data must be real");
647 x = (RVector) data_to_vector(s1, RE);
648 y = (RVector) data_to_vector(s2, RE);
649 z = (RMatrix) data_to_matrix(mat, RE);
651 xlsave1(result);
652 result = NIL;
653 for (i = 0; i < n - 1; i++) {
654 for (j = 0; j < m - 1; j++) {
655 result = add_contour_point(i, j, i, j+1, x, y, z, v, result);
656 result = add_contour_point(i, j+1, i+1, j+1, x, y, z, v, result);
657 result = add_contour_point(i+1, j+1, i+1, j, x, y, z, v, result);
658 result = add_contour_point(i+1, j, i, j, x, y, z, v, result);
661 xlpop();
663 free_vector(x);
664 free_vector(y);
665 free_matrix(z, n);
667 return(result);
672 ;;; FFT
674 ;;; FIXME:ajr
675 ;;; ??replace with matlisp:fft and matlisp:ifft (the latter for inverse mapping)
677 (defun fft (x &optional inverse)
678 "Args: (x &optional inverse)
679 Returns unnormalized Fourier transform of X, or inverse transform if INVERSE
680 is true."
681 (check-sequence x)
682 (let* ((n (length x))
683 (mode (la-data-mode x))
684 (isign (if inverse -1 1))
685 (result (if (consp x) (make-list n) (make-array n))))
686 (let ((px (la-data-to-vector x mode-cx))
687 (work (la-vector (+ (* 4 n) 15) mode-re)))
688 (unwind-protect
689 (progn
690 (fft-front n px work isign)
691 (la-vector-to-data px n mode-cx result))
692 (la-free-vector px)
693 (la-free-vector work))
694 result)))
697 ;;; SWEEP Operator: FIXME: use matlisp
700 (defun make-sweep-front (x y w n p mode has_w x_mean result)
701 (declare (fixnum n p mode has_w))
702 (let ((x_data nil)
703 (result_data nil)
704 (val 0.0)
705 (dxi 0.0)
706 (dyi 0.0)
707 (dv 0.0)
708 (dw 0.0)
709 (sum_w 0.0)
710 (dxik 0.0)
711 (dxjk 0.0)
712 (dyj 0.0)
713 (dx_meani 0.0)
714 (dx_meanj 0.0)
715 (dy_mean 0.0)
716 (has-w (if (/= 0 has_w) t nil))
717 (RE 1))
718 (declare (long-float val dxi dyi dv dw sum_w dxik dxjk dyj
719 dx_meani dx_meanj dy_mean))
720 ;; (declare-double val dxi dyi dv dw sum_w dxik dxjk dyj
721 ;; dx_meani dx_meanj dy_mean)
723 (if (> mode RE) (error "not supported for complex data yet"))
725 (setf x_data (compound-data-seq x))
726 (setf result_data (compound-data-seq result))
728 ;; find the mean of y
729 (setf val 0.0)
730 (setf sum_w 0.0)
731 (dotimes (i n)
732 (declare (fixnum i))
733 (setf dyi (makedouble (aref y i)))
734 (when has-w
735 (setf dw (makedouble (aref w i)))
736 (incf sum_w dw)
737 (setf dyi (* dyi dw)))
738 (incf val dyi))
739 (if (not has-w) (setf sum_w (float n 0.0)))
740 (if (<= sum_w 0.0) (error "non positive sum of weights"))
741 (setf dy_mean (/ val sum_w))
743 ;; find the column means
744 (dotimes (j p)
745 (declare (fixnum j))
746 (setf val 0.0)
747 (dotimes (i n)
748 (declare (fixnum i))
749 (setf dxi (makedouble (aref x_data (+ (* p i) j))))
750 (when has-w
751 (setf dw (makedouble (aref w i)))
752 (setf dxi (* dxi dw)))
753 (incf val dxi))
754 (setf (aref x_mean j) (/ val sum_w)))
756 ;; put 1/sum_w in topleft, means on left, minus means on top
757 (setf (aref result_data 0) (/ 1.0 sum_w))
758 (dotimes (i p)
759 (declare (fixnum i))
760 (setf dxi (makedouble (aref x_mean i)))
761 (setf (aref result_data (+ i 1)) (- dxi))
762 (setf (aref result_data (* (+ i 1) (+ p 2))) dxi))
763 (setf (aref result_data (+ p 1)) (- dy_mean))
764 (setf (aref result_data (* (+ p 1) (+ p 2))) dy_mean)
766 ;; put sums of adjusted cross products in body
767 (dotimes (i p)
768 (declare (fixnum i))
769 (dotimes (j p)
770 (declare (fixnum j))
771 (setf val 0.0)
772 (dotimes (k n)
773 (declare (fixnum k))
774 (setf dxik (makedouble (aref x_data (+ (* p k) i))))
775 (setf dxjk (makedouble (aref x_data (+ (* p k) j))))
776 (setf dx_meani (makedouble (aref x_mean i)))
777 (setf dx_meanj (makedouble (aref x_mean j)))
778 (setf dv (* (- dxik dx_meani) (- dxjk dx_meanj)))
779 (when has-w
780 (setf dw (makedouble (aref w k)))
781 (setf dv (* dv dw)))
782 (incf val dv))
783 (setf (aref result_data (+ (* (+ i 1) (+ p 2)) (+ j 1))) val)
784 (setf (aref result_data (+ (* (+ j 1) (+ p 2)) (+ i 1))) val))
785 (setf val 0.0)
786 (dotimes (j n)
787 (declare (fixnum j))
788 (setf dxik (makedouble (aref x_data (+ (* p j) i))))
789 (setf dyj (makedouble (aref y j)))
790 (setf dx_meani (makedouble (aref x_mean i)))
791 (setf dv (* (- dxik dx_meani) (- dyj dy_mean)))
792 (when has-w
793 (setf dw (makedouble (aref w j)))
794 (setf dv (* dv dw)))
795 (incf val dv))
796 (setf (aref result_data (+ (* (+ i 1) (+ p 2)) (+ p 1))) val)
797 (setf (aref result_data (+ (* (+ p 1) (+ p 2)) (+ i 1))) val))
798 (setf val 0.0)
799 (dotimes (j n)
800 (declare (fixnum j))
801 (setf dyj (makedouble (aref y j)))
802 (setf dv (* (- dyj dy_mean) (- dyj dy_mean)))
803 (when has-w
804 (setf dw (makedouble (aref w j)))
805 (setf dv (* dv dw)))
806 (incf val dv))
807 (setf (aref result_data (+ (* (+ p 1) (+ p 2)) (+ p 1))) val)))
809 ;;; FIXME: use matlisp
810 (defun sweep-in-place-front (a rows cols mode k tol)
811 "Sweep algorithm for linear regression."
812 (declare (long-float tol))
813 (declare (fixnum rows cols mode k))
814 (let ((data nil)
815 (pivot 0.0)
816 (aij 0.0)
817 (aik 0.0)
818 (akj 0.0)
819 (akk 0.0)
820 (RE 1))
821 (declare (long-float pivot aij aik akj akk))
823 (if (> mode RE) (error "not supported for complex data yet"))
824 (if (or (< k 0) (>= k rows) (>= k cols)) (error "index out of range"))
826 (setf tol (max tol machine-epsilon))
827 (setf data (compound-data-seq a))
829 (setf pivot (makedouble (aref data (+ (* cols k) k))))
831 (cond
832 ((or (> pivot tol) (< pivot (- tol)))
833 (dotimes (i rows)
834 (declare (fixnum i))
835 (dotimes (j cols)
836 (declare (fixnum j))
837 (when (and (/= i k) (/= j k))
838 (setf aij (makedouble (aref data (+ (* cols i) j))))
839 (setf aik (makedouble (aref data (+ (* cols i) k))))
840 (setf akj (makedouble (aref data (+ (* cols k) j))))
841 (setf aij (- aij (/ (* aik akj) pivot)))
842 (setf (aref data (+ (* cols i) j)) aij))))
844 (dotimes (i rows)
845 (declare (fixnum i))
846 (setf aik (makedouble (aref data (+ (* cols i) k))))
847 (when (/= i k)
848 (setf aik (/ aik pivot))
849 (setf (aref data (+ (* cols i) k)) aik)))
851 (dotimes (j cols)
852 (declare (fixnum j))
853 (setf akj (makedouble (aref data (+ (* cols k) j))))
854 (when (/= j k)
855 (setf akj (- (/ akj pivot)))
856 (setf (aref data (+ (* cols k) j)) akj)))
858 (setf akk (/ 1.0 pivot))
859 (setf (aref data (+ (* cols k) k)) akk)
861 (t 0))))
863 ;; FIXME: use matlisp
864 (defun make-sweep-matrix (x y &optional w)
865 "Args: (x y &optional weights)
866 X is matrix, Y and WEIGHTS are sequences. Returns the sweep matrix of the
867 (weighted) regression of Y on X"
868 (check-matrix x)
869 (check-sequence y)
870 (if w (check-sequence w))
871 (let ((n (num-rows x))
872 (p (num-cols x)))
873 (if (/= n (length y)) (error "dimensions do not match"))
874 (if (and w (/= n (length w))) (error "dimensions do not match"))
875 (let ((mode (max (la-data-mode x)
876 (la-data-mode x)
877 (if w (la-data-mode w) 0)))
878 (result (make-array (list (+ p 2) (+ p 2))))
879 (x-mean (make-array p))
880 (y (coerce y 'vector))
881 (w (if w (coerce w 'vector)))
882 (has-w (if w 1 0)))
883 (make-sweep-front x y w n p mode has-w x-mean result)
884 result)))
886 (defun sweep-in-place (a k tol)
887 (check-matrix a)
888 (check-one-fixnum k)
889 (check-one-real tol)
890 (let ((rows (num-rows a))
891 (cols (num-cols a))
892 (mode (la-data-mode a)))
893 (let ((swept (sweep-in-place-front a rows cols mode k tol)))
894 (if (/= 0 swept) t nil))))
896 (defun sweep-operator (a columns &optional tolerances)
897 "Args: (a indices &optional tolerances)
898 A is a matrix, INDICES a sequence of the column indices to be swept. Returns
899 a list of the swept result and the list of the columns actually swept. (See
900 MULTREG documentation.) If supplied, TOLERANCES should be a list of real
901 numbers the same length as INDICES. An index will only be swept if its pivot
902 element is larger than the corresponding element of TOLERANCES."
903 (check-matrix a)
904 (check-sequence columns)
905 (if tolerances (check-sequence tolerances))
906 (check-real a)
907 (check-fixnum columns)
908 (if tolerances (check-real tolerances))
909 (do ((tol .0000001)
910 (result (copy-array a))
911 (swept-columns nil)
912 (columns (coerce columns 'list) (cdr columns))
913 (tolerances (if (consp tolerances) (coerce tolerances 'list))
914 (if (consp tolerances) (cdr tolerances))))
915 ((null columns) (list result swept-columns))
916 (let ((col (first columns))
917 (tol (if (consp tolerances) (first tolerances) tol)))
918 (if (sweep-in-place result col tol)
919 (setf swept-columns (cons col swept-columns))))))
923 ;;; AX+Y
926 ;;; matlisp:axpy
928 (defun ax+y (a x y &optional lower)
929 "Args (a x y &optional lower)
930 Returns (+ (matmult A X) Y). If LOWER is not nil, A is taken to be lower
931 triangular.
932 This could probably be made more efficient."
933 (check-square-matrix a)
934 (check-sequence x)
935 (check-sequence y)
936 (check-real a)
937 (check-real x)
938 (check-real y)
939 (let* ((n (num-rows a))
940 (result (make-list n))
941 (a (compound-data-seq a)))
942 (declare (fixnum n))
943 (if (or (/= n (length x)) (/= n (length y)))
944 (error "dimensions do not match"))
945 (do* ((tx (make-next-element x) (make-next-element x))
946 (ty (make-next-element y))
947 (tr (make-next-element result))
948 (i 0 (+ i 1))
949 (start 0 (+ start n))
950 (end (if lower (+ i 1) n) (if lower (+ i 1) n)))
951 ((<= n i) result)
952 (declare (fixnum i start end))
953 (let ((val (get-next-element ty i)))
954 (dotimes (j end)
955 (declare (fixnum j))
956 (setf val (+ val (* (get-next-element tx j)
957 (aref a (+ start j))))))
958 (set-next-element tr i val)))))