lsfloat needed for mach-eps.
[CommonLispStat.git] / linalg.lsp
blobdc9874117ee50c68aa00fff7a6dfda672d2a3300
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 (in-package :cl-user)
30 (defpackage :lisp-stat-linalg
31 (:use :common-lisp
32 :cffi
33 :lisp-stat-ffi-int
34 :lisp-stat-math
35 :lisp-stat-types
36 :lisp-stat-float
37 :lisp-stat-compound-data
38 :lisp-stat-linalg-data
39 :lisp-stat-matrix)
40 (:shadowing-import-from :lisp-stat-math
41 expt + - * / ** mod rem abs 1+ 1- log exp sqrt sin cos tan
42 asin acos atan sinh cosh tanh asinh acosh atanh float random
43 truncate floor ceiling round minusp zerop plusp evenp oddp
44 < <= = /= >= > complex conjugate realpart imagpart phase
45 min max logand logior logxor lognot ffloor fceiling
46 ftruncate fround signum cis)
47 (:export chol-decomp lu-decomp lu-solve determinant inverse sv-decomp
48 qr-decomp rcondest make-rotation spline kernel-dens kernel-smooth
49 fft make-sweep-matrix sweep-operator ax+y eigen
51 check-real ;; for optimize
53 ;; the following are more matrix-centric (init from lsbasics, stats)
54 covariance-matrix matrix print-matrix solve
55 backsolve eigenvalues eigenvectors accumulate cumsum combine
56 lowess
60 (in-package #:lisp-stat-linalg)
62 ;;; CFFI Support
64 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
65 ;;;
66 ;;; Lisp Interfaces to Linear Algebra Routines
67 ;;;
68 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
70 ;;;
71 ;;; Cholesky Decomposition
72 ;;;
74 (cffi:defcfun ("ccl_chol_decomp_front" ccl-chol-decomp-front)
75 :int (x :pointer) (y :int) (z :pointer))
76 (defun chol-decomp-front (x y z)
77 (ccl-chol-decomp-front x y z))
79 ;;;;
80 ;;;; LU Decomposition
81 ;;;;
83 (cffi:defcfun ("ccl_lu_decomp_front" ccl-lu-decomp-front)
84 :int (x :pointer) (y :int) (z :pointer) (u :int) (v :pointer))
85 (defun lu-decomp-front (x y z u v)
86 (ccl-lu-decomp-front x y z u v))
88 (cffi:defcfun ("ccl_lu_solve_front" ccl-lu-solve-front)
89 :int (x :pointer) (y :int) (z :pointer) (u :pointer) (v :int))
90 (defun lu-solve-front (x y z u v)
91 (ccl-lu-solve-front x y z u v))
93 (cffi:defcfun ("ccl_lu_inverse_front" ccl-lu-inverse-front)
94 :int (x :pointer) (y :int) (z :pointer) (u :pointer) (v :int) (w :pointer))
95 (defun lu-inverse-front (x y z u v w)
96 (ccl-lu-inverse-front x y z u v w))
98 ;;;;
99 ;;;; SV Decomposition
100 ;;;;
102 (cffi:defcfun ("ccl_sv_decomp_front" ccl-sv-decomp-front)
103 :int (x :pointer) (y :int) (z :int) (u :pointer) (v :pointer))
104 (defun sv-decomp-front (x y z u v)
105 (ccl-sv-decomp-front x y z u v))
107 ;;;;
108 ;;;; QR Decomposition
109 ;;;;
111 (cffi:defcfun ("ccl_qr_decomp_front" ccl-qr-decomp-front)
112 :int (x :pointer) (y :int) (z :int) (u :pointer) (v :pointer) (w :int))
113 (defun qr-decomp-front (x y z u v w)
114 (ccl-qr-decomp-front x y z u v w))
116 ;;;;
117 ;;;; Estimate of Condition Number for Lower Triangular Matrix
118 ;;;;
120 (cffi:defcfun ("ccl_rcondest_front" ccl-rcondest-front)
121 :double (x :pointer) (y :int))
122 (defun rcondest-front (x y)
123 (ccl-rcondest-front x y))
125 ;;;;
126 ;;;; Make Rotation Matrix
127 ;;;;
129 (cffi:defcfun ("ccl_make_rotation_front" ccl-make-rotation-front)
130 :int (x :int) (y :pointer) (z :pointer) (u :pointer) (v :int) (w :double))
131 (defun make-rotation-front (x y z u v w)
132 (ccl-make-rotation-front x y z u v (float w 1d0)))
134 ;;;;
135 ;;;; Eigenvalues and Eigenvectors
136 ;;;;
138 (cffi:defcfun ("ccl_eigen_front" ccl-eigen-front)
139 :int (x :pointer) (y :int) (z :pointer) (u :pointer) (v :pointer))
140 (defun eigen-front (x y z u v)
141 (ccl-eigen-front x y z u v))
143 ;;;;
144 ;;;; Spline Interpolation
145 ;;;;
147 (cffi:defcfun ("ccl_range_to_rseq" ccl-range-to-rseq)
148 :int (x :int) (y :pointer) (z :int) (u :pointer))
149 (defun la-range-to-rseq (x y z u)
150 (ccl-range-to-rseq x y z u))
152 (cffi:defcfun ("ccl_spline_front" ccl-spline-front)
153 :int (x :int) (y :pointer) (z :pointer) (u :int) (v :pointer) (w :pointer) (a :pointer))
154 (defun spline-front (x y z u v w a)
155 (ccl-spline-front x y z u v w a))
157 ;;;;
158 ;;;; Kernel Density Estimators and Smoothers
159 ;;;;
161 (cffi:defcfun ("ccl_kernel_dens_front" ccl-kernel-dens-front)
162 :int (x :pointer) (y :int) (z :double) (u :pointer) (v :pointer) (w :int) (a :int))
163 (defun kernel-dens-front (x y z u v w a)
164 (ccl-kernel-dens-front x y (float z 1d0) u v w a))
166 (cffi:defcfun ("ccl_kernel_smooth_front" ccl-kernel-smooth-front)
167 :int (x :pointer) (y :pointer) (z :int) (u :double) (v :pointer) (w :pointer) (a :int) (b :int))
168 (defun kernel-smooth-front (x y z u v w a b)
169 (ccl-kernel-smooth-front x y z (float u 1d0) v w a b))
171 ;;;;
172 ;;;; Lowess Smoother Interface
173 ;;;;
175 (cffi:defcfun ("ccl_base_lowess_front" ccl-base-lowess-front)
176 :int (x :pointer) (y :pointer) (z :int) (u :double) (v :int) (w :double) (a :pointer) (b :pointer) (c :pointer))
177 (defun base-lowess-front (x y z u v w a b c)
178 (ccl-base-lowess-front x y z (float u 1d0) v (float w 1d0) a b c))
180 ;;;;
181 ;;;; FFT
182 ;;;;
184 (cffi:defcfun ("ccl_fft_front" ccl-fft-front)
185 :int (x :int) (y :pointer) (z :pointer) (u :int))
186 (defun fft-front (x y z u)
187 (ccl-fft-front x y z u))
192 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
193 ;;;;
194 ;;;; Lisp to C number conversion and checking
195 ;;;;
196 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
198 ;;;;
199 ;;;; Lisp to/from C sequence and matrix conversion and checking
200 ;;;;
202 (defun is-cons (a)
203 "FIXME:AJR this is not used anywhere?"
204 (if (consp a) 1 0))
206 (defun check-fixnum (a)
207 (if (/= 0 (la-data-mode a)) (error "not an integer sequence - ~s" a)))
209 (defun check-real (data)
210 (let ((data (compound-data-seq data)))
211 (cond
212 ((vectorp data)
213 (let ((n (length data)))
214 (declare (fixnum n))
215 (dotimes (i n)
216 (declare (fixnum i))
217 (check-one-real (aref data i)))))
218 ((consp data) (dolist (x data) (check-one-real x)))
219 (t (error "bad sequence - ~s" data)))))
221 (defun vec-assign (a i x) (setf (aref a i) x))
223 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
224 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
225 ;;;;
226 ;;;; Lisp Interfaces to Linear Algebra Routines
227 ;;;;
228 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
229 ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
231 ;;; FIXME: use dpbt[f2|rf], dpbstf, dpot[f2|rf]; dpptrf, zpbstf, zpbt[f2|rf]
232 ;;; remember: factorization = decomposition, depending on training.
234 (defun chol-decomp (a &optional (maxoffl 0.0))
235 "Args: (a)
236 Modified Cholesky decomposition. A should be a square, symmetric matrix.
237 Computes lower triangular matrix L such that L L^T = A + D where D is a
238 diagonal matrix. If A is strictly positive definite D will be zero.
239 Otherwise, D is as small as possible to make A + D numerically strictly
240 positive definite. Returns a list (L (max D))."
241 (check-square-matrix a)
242 (check-real a)
243 (let* ((n (array-dimension a 0))
244 (result (make-array (list n n)))
245 (dpars (list maxoffl 0.0)))
246 (check-real dpars)
247 (let ((mat (la-data-to-matrix a +mode-re+))
248 (dp (la-data-to-vector dpars +mode-re+)))
249 (unwind-protect
250 (progn
251 (chol-decomp-front mat n dp)
252 (la-matrix-to-data mat n n +mode-re+ result)
253 (la-vector-to-data dp 2 +mode-re+ dpars))
254 (la-free-matrix mat n)
255 (la-free-vector dp)))
256 (list result (second dpars))))
259 ;;; REPLACE with
260 ;;; (matlisp:lu M)
261 ;;; i.e. result use by:
262 ;;; (setf (values (lu-out1 lu-out2 lu-out3)) (matlisp:lu my-matrix))
263 ;;; for solution, ...
264 ;;; for lu-solve:
265 ;;; (matlisp:gesv a b &opt ipivot)
267 (defun lu-decomp (a)
268 "Args: (a)
269 A is a square matrix of numbers (real or complex). Computes the LU
270 decomposition of A and returns a list of the form (LU IV D FLAG), where
271 LU is a matrix with the L part in the lower triangle, the U part in the
272 upper triangle (the diagonal entries of L are taken to be 1), IV is a vector
273 describing the row permutation used, D is 1 if the number of permutations
274 is odd, -1 if even, and FLAG is T if A is numerically singular, NIL otherwise.
275 Used bu LU-SOLVE."
276 (check-square-matrix a)
277 (let* ((n (array-dimension a 0))
278 (mode (max +mode-re+ (la-data-mode a)))
279 (result (list (make-array (list n n)) (make-array n) nil nil)))
280 (let ((mat (la-data-to-matrix a mode))
281 (iv (la-vector n +mode-in+))
282 (d (la-vector 1 +mode-re+))
283 (singular 0))
284 (unwind-protect
285 (progn
286 (setf singular (lu-decomp-front mat n iv mode d))
287 (la-matrix-to-data mat n n mode (first result))
288 (la-vector-to-data iv n +mode-in+ (second result))
289 (setf (third result) (la-get-double d 0))
290 (setf (fourth result) (if (= singular 0.0) nil t)))
291 (la-free-matrix mat n)
292 (la-free-vector iv)
293 (la-free-vector d)))
294 result))
296 (defun lu-solve (lu lb)
297 "Args: (lu b)
298 LU is the result of (LU-DECOMP A) for a square matrix A, B is a sequence.
299 Returns the solution to the equation Ax = B. Signals an error if A is
300 singular."
301 (let ((la (first lu))
302 (lidx (second lu)))
303 (check-square-matrix la)
304 (check-sequence lidx)
305 (check-sequence lb)
306 (check-fixnum lidx)
307 (let* ((n (num-rows la))
308 (result (make-sequence (if (consp lb) 'list 'vector) n))
309 (a-mode (la-data-mode la))
310 (b-mode (la-data-mode lb)))
311 (if (/= n (length lidx)) (error "index sequence is wrong length"))
312 (if (/= n (length lb)) (error "right hand side is wrong length"))
313 (let* ((mode (max +mode-re+ a-mode b-mode))
314 (a (la-data-to-matrix la mode))
315 (indx (la-data-to-vector lidx +mode-in+))
316 (b (la-data-to-vector lb mode))
317 (singular 0))
318 (unwind-protect
319 (progn
320 (setf singular (lu-solve-front a n indx b mode))
321 (la-vector-to-data b n mode result))
322 (la-free-matrix a n)
323 (la-free-vector indx)
324 (la-free-vector b))
325 (if (/= 0.0 singular) (error "matrix is (numerically) singular"))
326 result))))
328 (defun determinant (a)
329 "Args: (m)
330 Returns the determinant of the square matrix M."
331 (let* ((lu (lu-decomp a))
332 (la (first lu))
333 (n (num-rows a))
334 (d1 (third lu))
335 (d2 0.d0))
336 (declare (fixnum n))
337 (flet ((fabs (x) (float (abs x) 0.d0)))
338 (dotimes (i n (* d1 (exp d2)))
339 (declare (fixnum i))
340 (let* ((x (aref la i i))
341 (magn (fabs x)))
342 (if (= 0.0 magn) (return 0.d0))
343 (setf d1 (* d1 (/ x magn)))
344 (setf d2 (+ d2 (log magn))))))))
346 (defun inverse (a)
347 "Args: (m)
348 Returns the inverse of the the square matrix M; signals an error if M is ill
349 conditioned or singular"
350 (check-square-matrix a)
351 (let ((n (num-rows a))
352 (mode (max +mode-re+ (la-data-mode a))))
353 (declare (fixnum n))
354 (let ((result (make-array (list n n) :initial-element 0)))
355 (dotimes (i n)
356 (declare (fixnum i))
357 (setf (aref result i i) 1))
358 (let ((mat (la-data-to-matrix a mode))
359 (inv (la-data-to-matrix result mode))
360 (iv (la-vector n +mode-in+))
361 (v (la-vector n mode))
362 (singular 0))
363 (unwind-protect
364 (progn
365 (setf singular (lu-inverse-front mat n iv v mode inv))
366 (la-matrix-to-data inv n n mode result))
367 (la-free-matrix mat n)
368 (la-free-matrix inv n)
369 (la-free-vector iv)
370 (la-free-vector v))
371 (if (/= singular 0) (error "matrix is (numerically) singular"))
372 result))))
374 ;;;;
375 ;;;; SV Decomposition
376 ;;;;
378 (defun sv-decomp (a)
379 "Args: (a)
380 A is a matrix of real numbers with at least as many rows as columns.
381 Computes the singular value decomposition of A and returns a list of the form
382 (U W V FLAG) where U and V are matrices whose columns are the left and right
383 singular vectors of A and W is the sequence of singular values of A. FLAG is T
384 if the algorithm converged, NIL otherwise."
385 (check-matrix a)
386 (let* ((m (num-rows a))
387 (n (num-cols a))
388 (mode (max +mode-re+ (la-data-mode a)))
389 (result (list (make-array (list m n))
390 (make-array n)
391 (make-array (list n n))
392 nil)))
393 (if (< m n) (error "number of rows less than number of columns"))
394 (if (= mode +mode-cx+) (error "complex SVD not available yet"))
395 (let ((mat (la-data-to-matrix a mode))
396 (w (la-vector n +mode-re+))
397 (v (la-matrix n n +mode-re+))
398 (converged 0))
399 (unwind-protect
400 (progn
401 (setf converged (sv-decomp-front mat m n w v))
402 (la-matrix-to-data mat m n mode (first result))
403 (la-vector-to-data w n mode (second result))
404 (la-matrix-to-data v n n mode (third result))
405 (setf (fourth result) (if (/= 0.0 converged) t nil)))
406 (la-free-matrix mat m)
407 (la-free-vector w)
408 (la-free-matrix v n))
409 result)))
412 ;;;;
413 ;;;; QR Decomposition
414 ;;;;
416 (defun qr-decomp (a &optional pivot)
417 "Args: (a &optional pivot)
418 A is a matrix of real numbers with at least as many rows as columns. Computes
419 the QR factorization of A and returns the result in a list of the form (Q R).
420 If PIVOT is true the columns of X are first permuted to place insure the
421 absolute values of the diagonal elements of R are nonincreasing. In this case
422 the result includes a third element, a list of the indices of the columns in
423 the order in which they were used."
424 (check-matrix a)
425 (let* ((m (num-rows a))
426 (n (num-cols a))
427 (mode (max +mode-re+ (la-data-mode a)))
428 (p (if pivot 1 0))
429 (result (if pivot
430 (list (make-array (list m n))
431 (make-array (list n n))
432 (make-array n))
433 (list (make-array (list m n)) (make-array (list n n))))))
434 (if (< m n) (error "number of rows less than number of columns"))
435 (if (= mode +mode-cx+) (error "complex QR decomposition not available yet"))
436 (let ((mat (la-data-to-matrix a mode))
437 (v (la-matrix n n +mode-re+))
438 (jpvt (la-vector n +mode-in+)))
439 (unwind-protect
440 (progn
441 (qr-decomp-front mat m n v jpvt p)
442 (la-matrix-to-data mat m n mode (first result))
443 (la-matrix-to-data v n n mode (second result))
444 (if pivot (la-vector-to-data jpvt n +mode-in+ (third result))))
445 (la-free-matrix mat m)
446 (la-free-matrix v n)
447 (la-free-vector jpvt))
448 result)))
450 ;;;;
451 ;;;; Estimate of Condition Number for Lower Triangular Matrix
452 ;;;;
454 (defun rcondest (a)
455 "Args: (a)
456 Returns an estimate of the reciprocal of the L1 condition number of an upper
457 triangular matrix a."
458 (check-square-matrix a)
459 (let ((mode (max +mode-re+ (la-data-mode a)))
460 (n (num-rows a)))
461 (if (= mode +mode-cx+)
462 (error "complex condition estimate not available yet"))
463 (let ((mat (la-data-to-matrix a mode))
464 (est 0.0))
465 (unwind-protect
466 (setf est (rcondest-front mat n))
467 (la-free-matrix mat n))
468 est)))
470 ;;;;
471 ;;;; Make Rotation Matrix
472 ;;;;
474 (defun make-rotation (x y &optional alpha)
475 "Args: (x y &optional alpha)
476 Returns a rotation matrix for rotating from X to Y, or from X toward Y
477 by angle ALPHA, in radians. X and Y are sequences of the same length."
478 (check-sequence x)
479 (check-sequence y)
480 (if alpha (check-one-real alpha))
481 (let* ((n (length x))
482 (mode (max +mode-re+ (la-data-mode x) (la-data-mode y)))
483 (use-angle (if alpha 1 0))
484 (angle (if alpha (float alpha 0.0) 0.0))
485 (result (make-array (list n n))))
486 (if (/= n (length y)) (error "sequences not the same length"))
487 (if (= mode +mode-cx+) (error "complex data not supported yet"))
488 (let ((px (la-data-to-vector x +mode-re+))
489 (py (la-data-to-vector y +mode-re+))
490 (rot (la-matrix n n +mode-re+)))
491 (unwind-protect
492 (progn
493 (make-rotation-front n rot px py use-angle angle)
494 (la-matrix-to-data rot n n +mode-re+ result))
495 (la-free-vector px)
496 (la-free-vector py)
497 (la-free-matrix rot n))
498 result)))
500 ;;;;
501 ;;;; Eigenvalues and Vectors
502 ;;;;
504 (defun eigen (a)
505 "Args: (a)
506 Returns list of list of eigenvalues and list of eigenvectors of square,
507 symmetric matrix A. Third element of result is NIL if algorithm converges.
508 If the algorithm does not converge, the third element is an integer I.
509 In this case the eigenvalues 0, ..., I are not reliable."
510 (check-square-matrix a)
511 (let ((mode (max +mode-re+ (la-data-mode a)))
512 (n (num-rows a)))
513 (if (= mode +mode-cx+) (error "matrix must be real and symmetric"))
514 (let ((evals (make-array n))
515 (evecs (make-list (* n n)))
516 (pa (la-data-to-vector (compound-data-seq a) +mode-re+))
517 (w (la-vector n +mode-re+))
518 (z (la-vector (* n n) +mode-re+))
519 (fv1 (la-vector n +mode-re+))
520 (ierr 0))
521 (unwind-protect
522 (progn
523 (setf ierr (eigen-front pa n w z fv1))
524 (la-vector-to-data w n +mode-re+ evals)
525 (la-vector-to-data z (* n n) +mode-re+ evecs))
526 (la-free-vector pa)
527 (la-free-vector z)
528 (la-free-vector w)
529 (la-free-vector fv1))
530 (list (nreverse evals)
531 (nreverse (mapcar #'(lambda (x) (coerce x 'vector))
532 (split-list evecs n)))
533 (if (/= 0 ierr) (- n ierr))))))
535 ;;;;
536 ;;;; Spline Interpolation
537 ;;;;
539 (defun make-smoother-args (x y xvals)
540 (check-sequence x)
541 (check-real x)
542 (when y
543 (check-sequence y)
544 (check-real y))
545 (unless (integerp xvals)
546 (check-sequence xvals)
547 (check-real xvals))
548 (let* ((n (length x))
549 (ns (if (integerp xvals) xvals (length xvals)))
550 (result (list (make-list ns) (make-list ns))))
551 (if (and y (/= n (length y))) (error "sequences not the same length"))
552 (list x y n (if (integerp xvals) 0 1) ns xvals result)))
554 (defun get-smoother-result (args) (seventh args))
556 (defmacro with-smoother-data ((x y xvals is-reg) &rest body)
557 `(progn
558 (check-sequence ,x)
559 (check-real ,x)
560 (when ,is-reg
561 (check-sequence ,y)
562 (check-real ,y))
563 (unless (integerp ,xvals)
564 (check-sequence ,xvals)
565 (check-real ,xvals))
566 (let* ((supplied (not (integerp ,xvals)))
567 (n (length ,x))
568 (ns (if supplied (length ,xvals) ,xvals))
569 (result (list (make-list ns) (make-list ns))))
570 (if (and ,is-reg (/= n (length ,y)))
571 (error "sequences not the same length"))
572 (if (and (not supplied) (< ns 2))
573 (error "too few points for interpolation"))
574 (let* ((px (la-data-to-vector ,x +mode-re+))
575 (py (if ,is-reg (la-data-to-vector ,y +mode-re+)))
576 (pxs (if supplied
577 (la-data-to-vector ,xvals +mode-re+)
578 (la-vector ns +mode-re+)))
579 (pys (la-vector ns +mode-re+)))
580 (unless supplied (la-range-to-rseq n px ns pxs))
581 (unwind-protect
582 (progn ,@body
583 (la-vector-to-data pxs ns +mode-re+ (first result))
584 (la-vector-to-data pys ns +mode-re+ (second result)))
585 (la-free-vector px)
586 (if ,is-reg (la-free-vector py))
587 (la-free-vector pxs)
588 (la-free-vector pys))
589 result))))
591 (defun spline (x y &key (xvals 30))
592 "Args: (x y &key xvals)
593 Returns list of x and y values of natural cubic spline interpolation of (X,Y).
594 X must be strictly increasing. XVALS can be an integer, the number of equally
595 spaced points to use in the range of X, or it can be a sequence of points at
596 which to interpolate."
597 (with-smoother-data (x y xvals t)
598 (let ((work (la-vector (* 2 n) +mode-re+))
599 (error 0))
600 (unwind-protect
601 (setf error (spline-front n px py ns pxs pys work))
602 (la-free-vector work))
603 (if (/= error 0) (error "bad data for splines")))))
605 ;;;;
606 ;;;; Kernel Density Estimators and Smoothers
607 ;;;;
609 (defun kernel-type-code (type)
610 (cond ((eq type 'u) 0)
611 ((eq type 't) 1)
612 ((eq type 'g) 2)
613 (t 3)))
615 (defun kernel-dens (x &key (type 'b) (width -1.0) (xvals 30))
616 "Args: (x &key xvals width type)
617 Returns list of x and y values of kernel density estimate of X. XVALS can be an
618 integer, the number of equally spaced points to use in the range of X, or it
619 can be a sequence of points at which to interpolate. WIDTH specifies the
620 window width. TYPE specifies the lernel and should be one of the symbols G, T,
621 U or B for gaussian, triangular, uniform or bisquare. The default is B."
622 (check-one-real width)
623 (with-smoother-data (x nil xvals nil) ;; warning about deleting unreachable code is TRUE -- 2nd arg=nil!
624 (let ((code (kernel-type-code type))
625 (error 0))
626 (setf error (kernel-dens-front px n width pxs pys ns code))
627 (if (/= 0 error) (error "bad kernel density data")))))
629 (defun kernel-smooth (x y &key (type 'b) (width -1.0) (xvals 30))
630 "Args: (x y &key xvals width type)
631 Returns list of x and y values of kernel smooth of (X,Y). XVALS can be an
632 integer, the number of equally spaced points to use in the range of X, or it
633 can be a sequence of points at which to interpolate. WIDTH specifies the
634 window width. TYPE specifies the lernel and should be one of the symbols G, T,
635 U or B for Gaussian, triangular, uniform or bisquare. The default is B."
636 (check-one-real width)
637 (with-smoother-data (x y xvals t)
638 (let ((code (kernel-type-code type))
639 (error 0))
640 (kernel-smooth-front px py n width pxs pys ns code)
641 ;; if we get the Lisp version ported from C, uncomment below and
642 ;; comment above. (thanks to Carlos Ungil for the initial CFFI
643 ;; work).
644 ;;(kernel-smooth-Cport px py n width pxs pys ns code)
645 (if (/= 0 error) (error "bad kernel density data")))))
649 (defun kernel-smooth-Cport (px py n width ;;wts wds ;; see above for mismatch?
650 xs ys ns ktype)
651 "Port of kernel_smooth (Lib/kernel.c) to Lisp.
652 FIXME:kernel-smooth-Cport : This is broken.
653 Until this is fixed, we are using Luke's C code and CFFI as glue."
654 (cond ((< n 1) 1.0)
655 ((and (< n 2) (<= width 0)) 1.0)
656 (t (let* ((xmin (min px))
657 (xmax (max px))
658 (width (/ (- xmax xmin) (+ 1.0 (log n)))))
659 (dotimes (i (- ns 1))
660 (setf (aref ys i)
661 (let ((wsum 0.0)
662 (ysum 0.0))
663 (dotimes (j (- n 1)) )
664 ;;;possible nasty errors...
666 (let*
667 ((lwidth (if wds (* width (aref wds j)) width))
668 (lwt (* (kernel-Cport (aref xs i) (aref px j) lwidth ktype) ;; px?
669 (if wts (aref wts j) 1.0))))
670 (setf wsum (+ wsum lwt))
671 (setf ysum (if py (+ ysum (* lwt (aref py j)))))) ;; py? y?
673 ;;; end of errors
674 (if py
675 (if (> wsum 0.0)
676 (/ ysum wsum)
677 0.0)
678 (/ wsum n)))))
679 (values ys)))))
681 (defun kernel-Cport (x y w ktype)
682 "Port of kernel() (Lib/kernel.c) to Lisp.
683 x,y,w are doubles, type is an integer"
684 (if (<= w 0.0)
686 (let ((z (- x y)))
687 (cond ((eq ktype "B")
688 (let* ((w (* w 2.0))
689 (z (* z 0.5)))
690 (if (and (> z -0.5)
691 (< z 0.5))
692 (/ (/ (* 15.0 (* (- 1.0 (* 4 z z)) ;; k/w
693 (- 1.0 (* 4 z z)))) ;; k/w
694 8.0)
696 0)))
697 ((eq ktype "G")
698 (let* ((w (* w 0.25))
699 (z (* z 4.0))
700 (k (/ (exp (* -0.5 z z))
701 (sqrt (* 2 PI)))))
702 (/ k w)))
703 ((eq ktype "U")
704 (let* ((w (* 1.5 w))
705 (z (* z 0.75))
706 (k (if (< (abs z) 0.5)
708 0.0)))
709 (/ k w)))
710 ((eq ktype "T")
711 (cond ((and (> z -1.0)
712 (< z 0.0))
713 (+ 1.0 z)) ;; k
714 ((and (> z 0.0)
715 (< z 1.0))
716 (- 1.0 z)) ;; k
717 (t 0.0)))
718 (t (values 0.0))))))
721 ;;;;
722 ;;;; Lowess Smoother Interface
723 ;;;;
725 (defun |base-lowess| (s1 s2 f nsteps delta)
726 (check-sequence s1)
727 (check-sequence s2)
728 (check-real s1)
729 (check-real s2)
730 (check-one-real f)
731 (check-one-fixnum nsteps)
732 (check-one-real delta)
733 (let* ((n (length s1))
734 (result (make-list n)))
735 (if (/= n (length s2)) (error "sequences not the same length"))
736 (let ((x (la-data-to-vector s1 +mode-re+))
737 (y (la-data-to-vector s2 +mode-re+))
738 (ys (la-vector n +mode-re+))
739 (rw (la-vector n +mode-re+))
740 (res (la-vector n +mode-re+))
741 (error 0))
742 (unwind-protect
743 (progn
744 (setf error (base-lowess-front x y n f nsteps delta ys rw res))
745 (la-vector-to-data ys n +mode-re+ result))
746 (la-free-vector x)
747 (la-free-vector y)
748 (la-free-vector ys)
749 (la-free-vector rw)
750 (la-free-vector res))
751 (if (/= error 0) (error "bad data for lowess"))
752 result)))
755 static LVAL add_contour_point(i, j, k, l, x, y, z, v, result)
756 int i, j, k, l;
757 RVector x, y;
758 RMatrix z;
759 double v;
760 LVAL result;
762 LVAL pt;
763 double p, q;
765 if ((z[i][j] <= v && v < z[k][l]) || (z[k][l] <= v && v < z[i][j])) {
766 xlsave(pt);
767 pt = mklist(2, NIL);
768 p = (v - z[i][j]) / (z[k][l] - z[i][j]);
769 q = 1.0 - p;
770 rplaca(pt, cvflonum((FLOTYPE) (q * x[i] + p * x[k])));
771 rplaca(cdr(pt), cvflonum((FLOTYPE) (q * y[j] + p * y[l])));
772 result = cons(pt, result);
773 xlpop();
775 return(result);
778 LVAL xssurface_contour()
780 LVAL s1, s2, mat, result;
781 RVector x, y;
782 RMatrix z;
783 double v;
784 int i, j, n, m;
786 s1 = xsgetsequence();
787 s2 = xsgetsequence();
788 mat = xsgetmatrix();
789 v = makedouble(xlgetarg());
790 xllastarg();
792 n = seqlen(s1); m = seqlen(s2);
793 if (n != numrows(mat) || m != numcols(mat)) xlfail("dimensions do not match");
794 if (data_mode(s1) == CX || data_mode(s2) == CX || data_mode(mat) == CX)
795 xlfail("data must be real");
797 x = (RVector) data_to_vector(s1, RE);
798 y = (RVector) data_to_vector(s2, RE);
799 z = (RMatrix) data_to_matrix(mat, RE);
801 xlsave1(result);
802 result = NIL;
803 for (i = 0; i < n - 1; i++) {
804 for (j = 0; j < m - 1; j++) {
805 result = add_contour_point(i, j, i, j+1, x, y, z, v, result);
806 result = add_contour_point(i, j+1, i+1, j+1, x, y, z, v, result);
807 result = add_contour_point(i+1, j+1, i+1, j, x, y, z, v, result);
808 result = add_contour_point(i+1, j, i, j, x, y, z, v, result);
811 xlpop();
813 free_vector(x);
814 free_vector(y);
815 free_matrix(z, n);
817 return(result);
822 ;;; FFT
824 ;;; FIXME:ajr
825 ;;; ??replace with matlisp:fft and matlisp:ifft (the latter for inverse mapping)
827 (defun fft (x &optional inverse)
828 "Args: (x &optional inverse)
829 Returns unnormalized Fourier transform of X, or inverse transform if INVERSE
830 is true."
831 (check-sequence x)
832 (let* ((n (length x))
833 (mode (la-data-mode x))
834 (isign (if inverse -1 1))
835 (result (if (consp x) (make-list n) (make-array n))))
836 (let ((px (la-data-to-vector x +mode-cx+))
837 (work (la-vector (+ (* 4 n) 15) +mode-re+)))
838 (unwind-protect
839 (progn
840 (fft-front n px work isign)
841 (la-vector-to-data px n +mode-cx+ result))
842 (la-free-vector px)
843 (la-free-vector work))
844 result)))
847 ;;; SWEEP Operator: FIXME: use matlisp
850 (defun make-sweep-front (x y w n p mode has_w x_mean result)
851 (declare (fixnum n p mode has_w))
852 (let ((x_data nil)
853 (result_data nil)
854 (val 0.0)
855 (dxi 0.0)
856 (dyi 0.0)
857 (dv 0.0)
858 (dw 0.0)
859 (sum_w 0.0)
860 (dxik 0.0)
861 (dxjk 0.0)
862 (dyj 0.0)
863 (dx_meani 0.0)
864 (dx_meanj 0.0)
865 (dy_mean 0.0)
866 (has-w (if (/= 0 has_w) t nil))
867 (RE 1))
868 (declare (long-float val dxi dyi dv dw sum_w dxik dxjk dyj
869 dx_meani dx_meanj dy_mean))
870 ;; (declare-double val dxi dyi dv dw sum_w dxik dxjk dyj
871 ;; dx_meani dx_meanj dy_mean)
873 (if (> mode RE) (error "not supported for complex data yet"))
875 (setf x_data (compound-data-seq x))
876 (setf result_data (compound-data-seq result))
878 ;; find the mean of y
879 (setf val 0.0)
880 (setf sum_w 0.0)
881 (dotimes (i n)
882 (declare (fixnum i))
883 (setf dyi (makedouble (aref y i)))
884 (when has-w
885 (setf dw (makedouble (aref w i)))
886 (incf sum_w dw)
887 (setf dyi (* dyi dw)))
888 (incf val dyi))
889 (if (not has-w) (setf sum_w (float n 0.0)))
890 (if (<= sum_w 0.0) (error "non positive sum of weights"))
891 (setf dy_mean (/ val sum_w))
893 ;; find the column means
894 (dotimes (j p)
895 (declare (fixnum j))
896 (setf val 0.0)
897 (dotimes (i n)
898 (declare (fixnum i))
899 (setf dxi (makedouble (aref x_data (+ (* p i) j))))
900 (when has-w
901 (setf dw (makedouble (aref w i)))
902 (setf dxi (* dxi dw)))
903 (incf val dxi))
904 (setf (aref x_mean j) (/ val sum_w)))
906 ;; put 1/sum_w in topleft, means on left, minus means on top
907 (setf (aref result_data 0) (/ 1.0 sum_w))
908 (dotimes (i p)
909 (declare (fixnum i))
910 (setf dxi (makedouble (aref x_mean i)))
911 (setf (aref result_data (+ i 1)) (- dxi))
912 (setf (aref result_data (* (+ i 1) (+ p 2))) dxi))
913 (setf (aref result_data (+ p 1)) (- dy_mean))
914 (setf (aref result_data (* (+ p 1) (+ p 2))) dy_mean)
916 ;; put sums of adjusted cross products in body
917 (dotimes (i p)
918 (declare (fixnum i))
919 (dotimes (j p)
920 (declare (fixnum j))
921 (setf val 0.0)
922 (dotimes (k n)
923 (declare (fixnum k))
924 (setf dxik (makedouble (aref x_data (+ (* p k) i))))
925 (setf dxjk (makedouble (aref x_data (+ (* p k) j))))
926 (setf dx_meani (makedouble (aref x_mean i)))
927 (setf dx_meanj (makedouble (aref x_mean j)))
928 (setf dv (* (- dxik dx_meani) (- dxjk dx_meanj)))
929 (when has-w
930 (setf dw (makedouble (aref w k)))
931 (setf dv (* dv dw)))
932 (incf val dv))
933 (setf (aref result_data (+ (* (+ i 1) (+ p 2)) (+ j 1))) val)
934 (setf (aref result_data (+ (* (+ j 1) (+ p 2)) (+ i 1))) val))
935 (setf val 0.0)
936 (dotimes (j n)
937 (declare (fixnum j))
938 (setf dxik (makedouble (aref x_data (+ (* p j) i))))
939 (setf dyj (makedouble (aref y j)))
940 (setf dx_meani (makedouble (aref x_mean i)))
941 (setf dv (* (- dxik dx_meani) (- dyj dy_mean)))
942 (when has-w
943 (setf dw (makedouble (aref w j)))
944 (setf dv (* dv dw)))
945 (incf val dv))
946 (setf (aref result_data (+ (* (+ i 1) (+ p 2)) (+ p 1))) val)
947 (setf (aref result_data (+ (* (+ p 1) (+ p 2)) (+ i 1))) val))
948 (setf val 0.0)
949 (dotimes (j n)
950 (declare (fixnum j))
951 (setf dyj (makedouble (aref y j)))
952 (setf dv (* (- dyj dy_mean) (- dyj dy_mean)))
953 (when has-w
954 (setf dw (makedouble (aref w j)))
955 (setf dv (* dv dw)))
956 (incf val dv))
957 (setf (aref result_data (+ (* (+ p 1) (+ p 2)) (+ p 1))) val)))
959 ;;; FIXME: use matlisp
960 (defun sweep-in-place-front (a rows cols mode k tol)
961 "Sweep algorithm for linear regression."
962 (declare (long-float tol))
963 (declare (fixnum rows cols mode k))
964 (let ((data nil)
965 (pivot 0.0)
966 (aij 0.0)
967 (aik 0.0)
968 (akj 0.0)
969 (akk 0.0)
970 (RE 1))
971 (declare (long-float pivot aij aik akj akk))
973 (if (> mode RE) (error "not supported for complex data yet"))
974 (if (or (< k 0) (>= k rows) (>= k cols)) (error "index out of range"))
976 (setf tol (max tol machine-epsilon))
977 (setf data (compound-data-seq a))
979 (setf pivot (makedouble (aref data (+ (* cols k) k))))
981 (cond
982 ((or (> pivot tol) (< pivot (- tol)))
983 (dotimes (i rows)
984 (declare (fixnum i))
985 (dotimes (j cols)
986 (declare (fixnum j))
987 (when (and (/= i k) (/= j k))
988 (setf aij (makedouble (aref data (+ (* cols i) j))))
989 (setf aik (makedouble (aref data (+ (* cols i) k))))
990 (setf akj (makedouble (aref data (+ (* cols k) j))))
991 (setf aij (- aij (/ (* aik akj) pivot)))
992 (setf (aref data (+ (* cols i) j)) aij))))
994 (dotimes (i rows)
995 (declare (fixnum i))
996 (setf aik (makedouble (aref data (+ (* cols i) k))))
997 (when (/= i k)
998 (setf aik (/ aik pivot))
999 (setf (aref data (+ (* cols i) k)) aik)))
1001 (dotimes (j cols)
1002 (declare (fixnum j))
1003 (setf akj (makedouble (aref data (+ (* cols k) j))))
1004 (when (/= j k)
1005 (setf akj (- (/ akj pivot)))
1006 (setf (aref data (+ (* cols k) j)) akj)))
1008 (setf akk (/ 1.0 pivot))
1009 (setf (aref data (+ (* cols k) k)) akk)
1011 (t 0))))
1013 ;; FIXME: use matlisp
1014 (defun make-sweep-matrix (x y &optional w)
1015 "Args: (x y &optional weights)
1016 X is matrix, Y and WEIGHTS are sequences. Returns the sweep matrix of the
1017 (weighted) regression of Y on X"
1018 (check-matrix x)
1019 (check-sequence y)
1020 (if w (check-sequence w))
1021 (let ((n (num-rows x))
1022 (p (num-cols x)))
1023 (if (/= n (length y)) (error "dimensions do not match"))
1024 (if (and w (/= n (length w))) (error "dimensions do not match"))
1025 (let ((mode (max (la-data-mode x)
1026 (la-data-mode x)
1027 (if w (la-data-mode w) 0)))
1028 (result (make-array (list (+ p 2) (+ p 2))))
1029 (x-mean (make-array p))
1030 (y (coerce y 'vector))
1031 (w (if w (coerce w 'vector)))
1032 (has-w (if w 1 0)))
1033 (make-sweep-front x y w n p mode has-w x-mean result)
1034 result)))
1036 (defun sweep-in-place (a k tol)
1037 (check-matrix a)
1038 (check-one-fixnum k)
1039 (check-one-real tol)
1040 (let ((rows (num-rows a))
1041 (cols (num-cols a))
1042 (mode (la-data-mode a)))
1043 (let ((swept (sweep-in-place-front a rows cols mode k tol)))
1044 (if (/= 0 swept) t nil))))
1046 (defun sweep-operator (a columns &optional tolerances)
1047 "Args: (a indices &optional tolerances)
1048 A is a matrix, INDICES a sequence of the column indices to be swept. Returns
1049 a list of the swept result and the list of the columns actually swept. (See
1050 MULTREG documentation.) If supplied, TOLERANCES should be a list of real
1051 numbers the same length as INDICES. An index will only be swept if its pivot
1052 element is larger than the corresponding element of TOLERANCES."
1053 (check-matrix a)
1054 (check-sequence columns)
1055 (if tolerances (check-sequence tolerances))
1056 (check-real a)
1057 (check-fixnum columns)
1058 (if tolerances (check-real tolerances))
1059 (do ((tol .0000001)
1060 (result (copy-array a))
1061 (swept-columns nil)
1062 (columns (coerce columns 'list) (cdr columns))
1063 (tolerances (if (consp tolerances) (coerce tolerances 'list))
1064 (if (consp tolerances) (cdr tolerances))))
1065 ((null columns) (list result swept-columns))
1066 (let ((col (first columns))
1067 (tol (if (consp tolerances) (first tolerances) tol)))
1068 (if (sweep-in-place result col tol)
1069 (setf swept-columns (cons col swept-columns))))))
1073 ;;; AX+Y
1076 ;;; matlisp:axpy
1078 (defun ax+y (a x y &optional lower)
1079 "Args (a x y &optional lower)
1080 Returns (+ (matmult A X) Y). If LOWER is not nil, A is taken to be lower
1081 triangular.
1082 This could probably be made more efficient."
1083 (check-square-matrix a)
1084 (check-sequence x)
1085 (check-sequence y)
1086 (check-real a)
1087 (check-real x)
1088 (check-real y)
1089 (let* ((n (num-rows a))
1090 (result (make-list n))
1091 (a (compound-data-seq a)))
1092 (declare (fixnum n))
1093 (if (or (/= n (length x)) (/= n (length y)))
1094 (error "dimensions do not match"))
1095 (do* ((tx (make-next-element x) (make-next-element x))
1096 (ty (make-next-element y))
1097 (tr (make-next-element result))
1098 (i 0 (+ i 1))
1099 (start 0 (+ start n))
1100 (end (if lower (+ i 1) n) (if lower (+ i 1) n)))
1101 ((<= n i) result)
1102 (declare (fixnum i start end))
1103 (let ((val (get-next-element ty i)))
1104 (dotimes (j end)
1105 (declare (fixnum j))
1106 (setf val (+ val (* (get-next-element tx j)
1107 (aref a (+ start j))))))
1108 (set-next-element tr i val)))))
1114 ;;;;
1115 ;;;; Linear Algebra Functions
1116 ;;;;
1118 (defun matrix (dim data)
1119 "Args: (dim data)
1120 returns a matrix of dimensions DIM initialized using sequence DATA
1121 in row major order."
1122 (let ((dim (coerce dim 'list))
1123 (data (coerce data 'list)))
1124 (make-array dim :initial-contents (split-list data (nth 1 dim)))))
1126 (defun flatsize (x)
1127 (length x)) ;; FIXME: defined badly!!
1129 (defun print-matrix (a &optional (stream *standard-output*))
1130 "Args: (matrix &optional stream)
1131 Prints MATRIX to STREAM in a nice form that is still machine readable"
1132 (unless (matrixp a) (error "not a matrix - ~a" a))
1133 (let ((size (min 15 (max (map-elements #'flatsize a)))))
1134 (format stream "#2a(~%")
1135 (dolist (x (row-list a))
1136 (format stream " (")
1137 (let ((n (length x)))
1138 (dotimes (i n)
1139 (let ((y (aref x i)))
1140 (cond
1141 ((integerp y) (format stream "~vd" size y))
1142 ((floatp y) (format stream "~vg" size y))
1143 (t (format stream "~va" size y))))
1144 (if (< i (- n 1)) (format stream " "))))
1145 (format stream ")~%"))
1146 (format stream " )~%")
1147 nil))
1149 (defun solve (a b)
1150 "Args: (a b)
1151 Solves A x = B using LU decomposition and backsolving. B can be a sequence
1152 or a matrix."
1153 (let ((lu (lu-decomp a)))
1154 (if (matrixp b)
1155 (apply #'bind-columns
1156 (mapcar #'(lambda (x) (lu-solve lu x)) (column-list b)))
1157 (lu-solve lu b))))
1159 (defun backsolve (a b)
1160 "Args: (a b)
1161 Solves A x = B by backsolving, assuming A is upper triangular. B must be a
1162 sequence. For use with qr-decomp."
1163 (let* ((n (length b))
1164 (sol (make-array n)))
1165 (dotimes (i n)
1166 (let* ((k (- n i 1))
1167 (val (elt b k)))
1168 (dotimes (j i)
1169 (let ((l (- n j 1)))
1170 (setq val (- val (* (aref sol l) (aref a k l))))))
1171 (setf (aref sol k) (/ val (aref a k k)))))
1172 (if (listp b) (coerce sol 'list) sol)))
1174 (defun eigenvalues (a)
1175 "Args: (a)
1176 Returns list of eigenvalues of square, symmetric matrix A"
1177 (first (eigen a)))
1179 (defun eigenvectors (a)
1180 "Args: (a)
1181 Returns list of eigenvectors of square, symmetric matrix A"
1182 (second (eigen a)))
1184 (defun accumulate (f s)
1185 "Args: (f s)
1186 Accumulates elements of sequence S using binary function F.
1187 (accumulate #'+ x) returns the cumulative sum of x."
1188 (let* ((result (list (elt s 0)))
1189 (tail result))
1190 (flet ((acc (dummy x)
1191 (rplacd tail (list (funcall f (first tail) x)))
1192 (setf tail (cdr tail))))
1193 (reduce #'acc s))
1194 (if (vectorp s) (coerce result 'vector) result)))
1196 (defun cumsum (x)
1197 "Args: (x)
1198 Returns the cumulative sum of X."
1199 (accumulate #'+ x))
1201 (defun combine (&rest args)
1202 "Args (&rest args)
1203 Returns sequence of elements of all arguments."
1204 (copy-seq (element-seq args)))
1206 (defun lowess (x y &key (f .25) (steps 2) (delta -1) sorted)
1207 "Args: (x y &key (f .25) (steps 2) delta sorted)
1208 Returns (list X YS) with YS the LOWESS fit. F is the fraction of data used for
1209 each point, STEPS is the number of robust iterations. Fits for points within
1210 DELTA of each other are interpolated linearly. If the X values setting SORTED
1211 to T speeds up the computation."
1212 (let ((x (if sorted x (sort-data x)))
1213 (y (if sorted y (select y (order x))))
1214 (delta (if (> delta 0.0) delta (/ (- (max x) (min x)) 50))))
1215 (list x)));; (|base-lowess| x y f steps delta))))