aboutsummaryrefslogtreecommitdiff
path: root/math/curfit/doc/curfit.hlp
blob: 35950c0872223ac163f895f855413b0e2da7174d (plain) (blame)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
.help curfit Jul84 "Math Package"
.ih
NAME
curfit -- curve fitting package
.ih
SYNOPSIS

.nf
        cvinit	(cv, curve_type, order, xmin, xmax)
	cvzero	(cv)
       cvaccum	(cv, x, y, weight, wtflag)
      cvreject	(cv, x, y, weight)
       cvsolve	(cv, ier)
	 cvfit	(cv, x, y, weight, npts, wtflag, ier)
       cvrefit	(cv, x, y, weight, ier)
    y = cveval	(cv, x)
      cvvector	(cv, x, yfit, npts)
       cvcoeff	(cv, coeff, ncoeff)
      cverrors	(cv, y, weight, yfit, rms, errors)
	cvsave	(cv, fit)
       cvstati  (cv, parameter, ival)
       cvstatr  (cv, parameter, ival)
     cvrestore	(cv, fit)
         cvset	(cv, curve_type, xmin, xmax, coeff, ncoeff)
        cvfree	(cv)
.fi
.ih
DESCRIPTION
The curfit package provides a set of routines for fitting data to functions
linear in their coefficients using least squares techniques. The numerical
technique employed is the solution of the normal equations by the
Cholesky method.
.ih
NOTES
The fitting function curve_type is chosen at run time from the following
list.

.nf
    LEGENDRE	# Legendre polynomials
    CHEBYSHEV	# Chebyshev polynomials
    SPLINE3	# cubic spline with uniformly spaced break points
    SPLINE1	# linear spline with uniformly spaced break points
.fi


The CURFIT package performs a weighted fit.
The weighting options are WTS_USR, WTS_UNIFORM and WTS_SPACING.
The user must supply a weight array. In WTS_UNIFORM mode the curfit
routines set the weights to 1. In WTS_USER mode the user must supply an
array of weight values.
In WTS_SPACING mode
the weights are set to the difference between adjacent data points.
The data must be sorted in x in order to use the WTS_SPACING mode.
In WTS_UNIFORM mode the reduced chi-squared returned by CVERRORS 
is the variance of the fit and the errors in the coefficients are scaled
by the square root of this variance. Otherwise the weights are
interpreted as one over the variance of the data and the true reduced
chi-squared is returned.

The routines assume that all the x values  of interest lie in the region
xmin <= x <= xmax. Checking for out of bounds x values is the responsibility
of the calling program. The package routines assume that INDEF values
have been removed from the data set prior to entering the package
routines.

In order to make the package definitions available to the calling program
an include <curfit.h> statement must be included in the user program.
CVINIT must be called before each fit. CVFREE frees space used by the
CURFIT package.
.ih
EXAMPLES
.nf
Example 1: Fit curve to data, unifrom weighting

    include <math/curfit.h>

    ...

    call cvinit (cv, CHEBYSHEV, 4, 1., 512.)

    call cvfit (cv, x, y, weight, 512, WTS_UNIFORM, ier)
    if (ier != OK)
	call error (...)

    do i = 1, 512 {
	x = i
	call printf ("%g %g\n")
	    call pargr (x)
	    call pargr (cveval (cv, x))
    }

    call cvfree (cv)


Example 2: Fit curve using accumulate mode, weight based on spacing

    include <math/curfit.h>

    ...

    old_x = x
    do i = 1, 512 {
	x = real (i)
	if (y[i] != INDEF) {
	    call cvaccum (cv, x, y, weight, x - old_x, WTS_USER)
	    old_x = x
	}
    }

    call cvsolve (cv, ier)
    if (ier != OK)
	call error (...)

    ...

    call cvfree (cv)


Example 3: Fit and subtract smooth curve from image lines

    include <math/curfit.h>

    ...

    call cvinit (cv, CHEBYSHEV, order, 1., 512.)

    do line = 1, nlines {
	inpix = imgl2r (im, line)
	outpix = impl2r (im, line)
	if (line == 1)
	    call cvfit (cv, x, Memr[inpix], weight, 512, WTS_UNIFORM, ier)
	else
	    call cvrefit (cv, x, Memr[inpix], weight, ier)
	if (ier != OK)
	    ...
	call cvvector (cv, x, y, 512)
	call asubr (Memr[inpix], y, Memr[outpix], 512)
    }

    call cvfree (cv)


Example 4: Fit curve, save fit for later use by CVEVAL. LEN_FIT must be a least
	   order + 7 elements long.

    include <math/curfit.h>

    real  fit[LEN_FIT]

    ...
    call cvinit (cv, CHEBYSHEV, order, xmin, xmax)
    call cvfit (cv, x, y, w, npts, WTS_UNIFORM, ier)
    if (ier != OK)
	...
    call cvsave (cv, fit)
    call cvfree (cv)
    ...
    call cvrestore (cv, fit)
    do i = 1, npts
	yfit[i] = cveval (cv, x[i])
    call cvfree (cv)
    ...
.fi