Chris@87
|
1 """
|
Chris@87
|
2 Objects for dealing with Laguerre series.
|
Chris@87
|
3
|
Chris@87
|
4 This module provides a number of objects (mostly functions) useful for
|
Chris@87
|
5 dealing with Laguerre series, including a `Laguerre` class that
|
Chris@87
|
6 encapsulates the usual arithmetic operations. (General information
|
Chris@87
|
7 on how this module represents and works with such polynomials is in the
|
Chris@87
|
8 docstring for its "parent" sub-package, `numpy.polynomial`).
|
Chris@87
|
9
|
Chris@87
|
10 Constants
|
Chris@87
|
11 ---------
|
Chris@87
|
12 - `lagdomain` -- Laguerre series default domain, [-1,1].
|
Chris@87
|
13 - `lagzero` -- Laguerre series that evaluates identically to 0.
|
Chris@87
|
14 - `lagone` -- Laguerre series that evaluates identically to 1.
|
Chris@87
|
15 - `lagx` -- Laguerre series for the identity map, ``f(x) = x``.
|
Chris@87
|
16
|
Chris@87
|
17 Arithmetic
|
Chris@87
|
18 ----------
|
Chris@87
|
19 - `lagmulx` -- multiply a Laguerre series in ``P_i(x)`` by ``x``.
|
Chris@87
|
20 - `lagadd` -- add two Laguerre series.
|
Chris@87
|
21 - `lagsub` -- subtract one Laguerre series from another.
|
Chris@87
|
22 - `lagmul` -- multiply two Laguerre series.
|
Chris@87
|
23 - `lagdiv` -- divide one Laguerre series by another.
|
Chris@87
|
24 - `lagval` -- evaluate a Laguerre series at given points.
|
Chris@87
|
25 - `lagval2d` -- evaluate a 2D Laguerre series at given points.
|
Chris@87
|
26 - `lagval3d` -- evaluate a 3D Laguerre series at given points.
|
Chris@87
|
27 - `laggrid2d` -- evaluate a 2D Laguerre series on a Cartesian product.
|
Chris@87
|
28 - `laggrid3d` -- evaluate a 3D Laguerre series on a Cartesian product.
|
Chris@87
|
29
|
Chris@87
|
30 Calculus
|
Chris@87
|
31 --------
|
Chris@87
|
32 - `lagder` -- differentiate a Laguerre series.
|
Chris@87
|
33 - `lagint` -- integrate a Laguerre series.
|
Chris@87
|
34
|
Chris@87
|
35 Misc Functions
|
Chris@87
|
36 --------------
|
Chris@87
|
37 - `lagfromroots` -- create a Laguerre series with specified roots.
|
Chris@87
|
38 - `lagroots` -- find the roots of a Laguerre series.
|
Chris@87
|
39 - `lagvander` -- Vandermonde-like matrix for Laguerre polynomials.
|
Chris@87
|
40 - `lagvander2d` -- Vandermonde-like matrix for 2D power series.
|
Chris@87
|
41 - `lagvander3d` -- Vandermonde-like matrix for 3D power series.
|
Chris@87
|
42 - `laggauss` -- Gauss-Laguerre quadrature, points and weights.
|
Chris@87
|
43 - `lagweight` -- Laguerre weight function.
|
Chris@87
|
44 - `lagcompanion` -- symmetrized companion matrix in Laguerre form.
|
Chris@87
|
45 - `lagfit` -- least-squares fit returning a Laguerre series.
|
Chris@87
|
46 - `lagtrim` -- trim leading coefficients from a Laguerre series.
|
Chris@87
|
47 - `lagline` -- Laguerre series of given straight line.
|
Chris@87
|
48 - `lag2poly` -- convert a Laguerre series to a polynomial.
|
Chris@87
|
49 - `poly2lag` -- convert a polynomial to a Laguerre series.
|
Chris@87
|
50
|
Chris@87
|
51 Classes
|
Chris@87
|
52 -------
|
Chris@87
|
53 - `Laguerre` -- A Laguerre series class.
|
Chris@87
|
54
|
Chris@87
|
55 See also
|
Chris@87
|
56 --------
|
Chris@87
|
57 `numpy.polynomial`
|
Chris@87
|
58
|
Chris@87
|
59 """
|
Chris@87
|
60 from __future__ import division, absolute_import, print_function
|
Chris@87
|
61
|
Chris@87
|
62 import warnings
|
Chris@87
|
63 import numpy as np
|
Chris@87
|
64 import numpy.linalg as la
|
Chris@87
|
65
|
Chris@87
|
66 from . import polyutils as pu
|
Chris@87
|
67 from ._polybase import ABCPolyBase
|
Chris@87
|
68
|
Chris@87
|
69 __all__ = [
|
Chris@87
|
70 'lagzero', 'lagone', 'lagx', 'lagdomain', 'lagline', 'lagadd',
|
Chris@87
|
71 'lagsub', 'lagmulx', 'lagmul', 'lagdiv', 'lagpow', 'lagval', 'lagder',
|
Chris@87
|
72 'lagint', 'lag2poly', 'poly2lag', 'lagfromroots', 'lagvander',
|
Chris@87
|
73 'lagfit', 'lagtrim', 'lagroots', 'Laguerre', 'lagval2d', 'lagval3d',
|
Chris@87
|
74 'laggrid2d', 'laggrid3d', 'lagvander2d', 'lagvander3d', 'lagcompanion',
|
Chris@87
|
75 'laggauss', 'lagweight']
|
Chris@87
|
76
|
Chris@87
|
77 lagtrim = pu.trimcoef
|
Chris@87
|
78
|
Chris@87
|
79
|
Chris@87
|
80 def poly2lag(pol):
|
Chris@87
|
81 """
|
Chris@87
|
82 poly2lag(pol)
|
Chris@87
|
83
|
Chris@87
|
84 Convert a polynomial to a Laguerre series.
|
Chris@87
|
85
|
Chris@87
|
86 Convert an array representing the coefficients of a polynomial (relative
|
Chris@87
|
87 to the "standard" basis) ordered from lowest degree to highest, to an
|
Chris@87
|
88 array of the coefficients of the equivalent Laguerre series, ordered
|
Chris@87
|
89 from lowest to highest degree.
|
Chris@87
|
90
|
Chris@87
|
91 Parameters
|
Chris@87
|
92 ----------
|
Chris@87
|
93 pol : array_like
|
Chris@87
|
94 1-D array containing the polynomial coefficients
|
Chris@87
|
95
|
Chris@87
|
96 Returns
|
Chris@87
|
97 -------
|
Chris@87
|
98 c : ndarray
|
Chris@87
|
99 1-D array containing the coefficients of the equivalent Laguerre
|
Chris@87
|
100 series.
|
Chris@87
|
101
|
Chris@87
|
102 See Also
|
Chris@87
|
103 --------
|
Chris@87
|
104 lag2poly
|
Chris@87
|
105
|
Chris@87
|
106 Notes
|
Chris@87
|
107 -----
|
Chris@87
|
108 The easy way to do conversions between polynomial basis sets
|
Chris@87
|
109 is to use the convert method of a class instance.
|
Chris@87
|
110
|
Chris@87
|
111 Examples
|
Chris@87
|
112 --------
|
Chris@87
|
113 >>> from numpy.polynomial.laguerre import poly2lag
|
Chris@87
|
114 >>> poly2lag(np.arange(4))
|
Chris@87
|
115 array([ 23., -63., 58., -18.])
|
Chris@87
|
116
|
Chris@87
|
117 """
|
Chris@87
|
118 [pol] = pu.as_series([pol])
|
Chris@87
|
119 deg = len(pol) - 1
|
Chris@87
|
120 res = 0
|
Chris@87
|
121 for i in range(deg, -1, -1):
|
Chris@87
|
122 res = lagadd(lagmulx(res), pol[i])
|
Chris@87
|
123 return res
|
Chris@87
|
124
|
Chris@87
|
125
|
Chris@87
|
126 def lag2poly(c):
|
Chris@87
|
127 """
|
Chris@87
|
128 Convert a Laguerre series to a polynomial.
|
Chris@87
|
129
|
Chris@87
|
130 Convert an array representing the coefficients of a Laguerre series,
|
Chris@87
|
131 ordered from lowest degree to highest, to an array of the coefficients
|
Chris@87
|
132 of the equivalent polynomial (relative to the "standard" basis) ordered
|
Chris@87
|
133 from lowest to highest degree.
|
Chris@87
|
134
|
Chris@87
|
135 Parameters
|
Chris@87
|
136 ----------
|
Chris@87
|
137 c : array_like
|
Chris@87
|
138 1-D array containing the Laguerre series coefficients, ordered
|
Chris@87
|
139 from lowest order term to highest.
|
Chris@87
|
140
|
Chris@87
|
141 Returns
|
Chris@87
|
142 -------
|
Chris@87
|
143 pol : ndarray
|
Chris@87
|
144 1-D array containing the coefficients of the equivalent polynomial
|
Chris@87
|
145 (relative to the "standard" basis) ordered from lowest order term
|
Chris@87
|
146 to highest.
|
Chris@87
|
147
|
Chris@87
|
148 See Also
|
Chris@87
|
149 --------
|
Chris@87
|
150 poly2lag
|
Chris@87
|
151
|
Chris@87
|
152 Notes
|
Chris@87
|
153 -----
|
Chris@87
|
154 The easy way to do conversions between polynomial basis sets
|
Chris@87
|
155 is to use the convert method of a class instance.
|
Chris@87
|
156
|
Chris@87
|
157 Examples
|
Chris@87
|
158 --------
|
Chris@87
|
159 >>> from numpy.polynomial.laguerre import lag2poly
|
Chris@87
|
160 >>> lag2poly([ 23., -63., 58., -18.])
|
Chris@87
|
161 array([ 0., 1., 2., 3.])
|
Chris@87
|
162
|
Chris@87
|
163 """
|
Chris@87
|
164 from .polynomial import polyadd, polysub, polymulx
|
Chris@87
|
165
|
Chris@87
|
166 [c] = pu.as_series([c])
|
Chris@87
|
167 n = len(c)
|
Chris@87
|
168 if n == 1:
|
Chris@87
|
169 return c
|
Chris@87
|
170 else:
|
Chris@87
|
171 c0 = c[-2]
|
Chris@87
|
172 c1 = c[-1]
|
Chris@87
|
173 # i is the current degree of c1
|
Chris@87
|
174 for i in range(n - 1, 1, -1):
|
Chris@87
|
175 tmp = c0
|
Chris@87
|
176 c0 = polysub(c[i - 2], (c1*(i - 1))/i)
|
Chris@87
|
177 c1 = polyadd(tmp, polysub((2*i - 1)*c1, polymulx(c1))/i)
|
Chris@87
|
178 return polyadd(c0, polysub(c1, polymulx(c1)))
|
Chris@87
|
179
|
Chris@87
|
180 #
|
Chris@87
|
181 # These are constant arrays are of integer type so as to be compatible
|
Chris@87
|
182 # with the widest range of other types, such as Decimal.
|
Chris@87
|
183 #
|
Chris@87
|
184
|
Chris@87
|
185 # Laguerre
|
Chris@87
|
186 lagdomain = np.array([0, 1])
|
Chris@87
|
187
|
Chris@87
|
188 # Laguerre coefficients representing zero.
|
Chris@87
|
189 lagzero = np.array([0])
|
Chris@87
|
190
|
Chris@87
|
191 # Laguerre coefficients representing one.
|
Chris@87
|
192 lagone = np.array([1])
|
Chris@87
|
193
|
Chris@87
|
194 # Laguerre coefficients representing the identity x.
|
Chris@87
|
195 lagx = np.array([1, -1])
|
Chris@87
|
196
|
Chris@87
|
197
|
Chris@87
|
198 def lagline(off, scl):
|
Chris@87
|
199 """
|
Chris@87
|
200 Laguerre series whose graph is a straight line.
|
Chris@87
|
201
|
Chris@87
|
202
|
Chris@87
|
203
|
Chris@87
|
204 Parameters
|
Chris@87
|
205 ----------
|
Chris@87
|
206 off, scl : scalars
|
Chris@87
|
207 The specified line is given by ``off + scl*x``.
|
Chris@87
|
208
|
Chris@87
|
209 Returns
|
Chris@87
|
210 -------
|
Chris@87
|
211 y : ndarray
|
Chris@87
|
212 This module's representation of the Laguerre series for
|
Chris@87
|
213 ``off + scl*x``.
|
Chris@87
|
214
|
Chris@87
|
215 See Also
|
Chris@87
|
216 --------
|
Chris@87
|
217 polyline, chebline
|
Chris@87
|
218
|
Chris@87
|
219 Examples
|
Chris@87
|
220 --------
|
Chris@87
|
221 >>> from numpy.polynomial.laguerre import lagline, lagval
|
Chris@87
|
222 >>> lagval(0,lagline(3, 2))
|
Chris@87
|
223 3.0
|
Chris@87
|
224 >>> lagval(1,lagline(3, 2))
|
Chris@87
|
225 5.0
|
Chris@87
|
226
|
Chris@87
|
227 """
|
Chris@87
|
228 if scl != 0:
|
Chris@87
|
229 return np.array([off + scl, -scl])
|
Chris@87
|
230 else:
|
Chris@87
|
231 return np.array([off])
|
Chris@87
|
232
|
Chris@87
|
233
|
Chris@87
|
234 def lagfromroots(roots):
|
Chris@87
|
235 """
|
Chris@87
|
236 Generate a Laguerre series with given roots.
|
Chris@87
|
237
|
Chris@87
|
238 The function returns the coefficients of the polynomial
|
Chris@87
|
239
|
Chris@87
|
240 .. math:: p(x) = (x - r_0) * (x - r_1) * ... * (x - r_n),
|
Chris@87
|
241
|
Chris@87
|
242 in Laguerre form, where the `r_n` are the roots specified in `roots`.
|
Chris@87
|
243 If a zero has multiplicity n, then it must appear in `roots` n times.
|
Chris@87
|
244 For instance, if 2 is a root of multiplicity three and 3 is a root of
|
Chris@87
|
245 multiplicity 2, then `roots` looks something like [2, 2, 2, 3, 3]. The
|
Chris@87
|
246 roots can appear in any order.
|
Chris@87
|
247
|
Chris@87
|
248 If the returned coefficients are `c`, then
|
Chris@87
|
249
|
Chris@87
|
250 .. math:: p(x) = c_0 + c_1 * L_1(x) + ... + c_n * L_n(x)
|
Chris@87
|
251
|
Chris@87
|
252 The coefficient of the last term is not generally 1 for monic
|
Chris@87
|
253 polynomials in Laguerre form.
|
Chris@87
|
254
|
Chris@87
|
255 Parameters
|
Chris@87
|
256 ----------
|
Chris@87
|
257 roots : array_like
|
Chris@87
|
258 Sequence containing the roots.
|
Chris@87
|
259
|
Chris@87
|
260 Returns
|
Chris@87
|
261 -------
|
Chris@87
|
262 out : ndarray
|
Chris@87
|
263 1-D array of coefficients. If all roots are real then `out` is a
|
Chris@87
|
264 real array, if some of the roots are complex, then `out` is complex
|
Chris@87
|
265 even if all the coefficients in the result are real (see Examples
|
Chris@87
|
266 below).
|
Chris@87
|
267
|
Chris@87
|
268 See Also
|
Chris@87
|
269 --------
|
Chris@87
|
270 polyfromroots, legfromroots, chebfromroots, hermfromroots,
|
Chris@87
|
271 hermefromroots.
|
Chris@87
|
272
|
Chris@87
|
273 Examples
|
Chris@87
|
274 --------
|
Chris@87
|
275 >>> from numpy.polynomial.laguerre import lagfromroots, lagval
|
Chris@87
|
276 >>> coef = lagfromroots((-1, 0, 1))
|
Chris@87
|
277 >>> lagval((-1, 0, 1), coef)
|
Chris@87
|
278 array([ 0., 0., 0.])
|
Chris@87
|
279 >>> coef = lagfromroots((-1j, 1j))
|
Chris@87
|
280 >>> lagval((-1j, 1j), coef)
|
Chris@87
|
281 array([ 0.+0.j, 0.+0.j])
|
Chris@87
|
282
|
Chris@87
|
283 """
|
Chris@87
|
284 if len(roots) == 0:
|
Chris@87
|
285 return np.ones(1)
|
Chris@87
|
286 else:
|
Chris@87
|
287 [roots] = pu.as_series([roots], trim=False)
|
Chris@87
|
288 roots.sort()
|
Chris@87
|
289 p = [lagline(-r, 1) for r in roots]
|
Chris@87
|
290 n = len(p)
|
Chris@87
|
291 while n > 1:
|
Chris@87
|
292 m, r = divmod(n, 2)
|
Chris@87
|
293 tmp = [lagmul(p[i], p[i+m]) for i in range(m)]
|
Chris@87
|
294 if r:
|
Chris@87
|
295 tmp[0] = lagmul(tmp[0], p[-1])
|
Chris@87
|
296 p = tmp
|
Chris@87
|
297 n = m
|
Chris@87
|
298 return p[0]
|
Chris@87
|
299
|
Chris@87
|
300
|
Chris@87
|
301 def lagadd(c1, c2):
|
Chris@87
|
302 """
|
Chris@87
|
303 Add one Laguerre series to another.
|
Chris@87
|
304
|
Chris@87
|
305 Returns the sum of two Laguerre series `c1` + `c2`. The arguments
|
Chris@87
|
306 are sequences of coefficients ordered from lowest order term to
|
Chris@87
|
307 highest, i.e., [1,2,3] represents the series ``P_0 + 2*P_1 + 3*P_2``.
|
Chris@87
|
308
|
Chris@87
|
309 Parameters
|
Chris@87
|
310 ----------
|
Chris@87
|
311 c1, c2 : array_like
|
Chris@87
|
312 1-D arrays of Laguerre series coefficients ordered from low to
|
Chris@87
|
313 high.
|
Chris@87
|
314
|
Chris@87
|
315 Returns
|
Chris@87
|
316 -------
|
Chris@87
|
317 out : ndarray
|
Chris@87
|
318 Array representing the Laguerre series of their sum.
|
Chris@87
|
319
|
Chris@87
|
320 See Also
|
Chris@87
|
321 --------
|
Chris@87
|
322 lagsub, lagmul, lagdiv, lagpow
|
Chris@87
|
323
|
Chris@87
|
324 Notes
|
Chris@87
|
325 -----
|
Chris@87
|
326 Unlike multiplication, division, etc., the sum of two Laguerre series
|
Chris@87
|
327 is a Laguerre series (without having to "reproject" the result onto
|
Chris@87
|
328 the basis set) so addition, just like that of "standard" polynomials,
|
Chris@87
|
329 is simply "component-wise."
|
Chris@87
|
330
|
Chris@87
|
331 Examples
|
Chris@87
|
332 --------
|
Chris@87
|
333 >>> from numpy.polynomial.laguerre import lagadd
|
Chris@87
|
334 >>> lagadd([1, 2, 3], [1, 2, 3, 4])
|
Chris@87
|
335 array([ 2., 4., 6., 4.])
|
Chris@87
|
336
|
Chris@87
|
337
|
Chris@87
|
338 """
|
Chris@87
|
339 # c1, c2 are trimmed copies
|
Chris@87
|
340 [c1, c2] = pu.as_series([c1, c2])
|
Chris@87
|
341 if len(c1) > len(c2):
|
Chris@87
|
342 c1[:c2.size] += c2
|
Chris@87
|
343 ret = c1
|
Chris@87
|
344 else:
|
Chris@87
|
345 c2[:c1.size] += c1
|
Chris@87
|
346 ret = c2
|
Chris@87
|
347 return pu.trimseq(ret)
|
Chris@87
|
348
|
Chris@87
|
349
|
Chris@87
|
350 def lagsub(c1, c2):
|
Chris@87
|
351 """
|
Chris@87
|
352 Subtract one Laguerre series from another.
|
Chris@87
|
353
|
Chris@87
|
354 Returns the difference of two Laguerre series `c1` - `c2`. The
|
Chris@87
|
355 sequences of coefficients are from lowest order term to highest, i.e.,
|
Chris@87
|
356 [1,2,3] represents the series ``P_0 + 2*P_1 + 3*P_2``.
|
Chris@87
|
357
|
Chris@87
|
358 Parameters
|
Chris@87
|
359 ----------
|
Chris@87
|
360 c1, c2 : array_like
|
Chris@87
|
361 1-D arrays of Laguerre series coefficients ordered from low to
|
Chris@87
|
362 high.
|
Chris@87
|
363
|
Chris@87
|
364 Returns
|
Chris@87
|
365 -------
|
Chris@87
|
366 out : ndarray
|
Chris@87
|
367 Of Laguerre series coefficients representing their difference.
|
Chris@87
|
368
|
Chris@87
|
369 See Also
|
Chris@87
|
370 --------
|
Chris@87
|
371 lagadd, lagmul, lagdiv, lagpow
|
Chris@87
|
372
|
Chris@87
|
373 Notes
|
Chris@87
|
374 -----
|
Chris@87
|
375 Unlike multiplication, division, etc., the difference of two Laguerre
|
Chris@87
|
376 series is a Laguerre series (without having to "reproject" the result
|
Chris@87
|
377 onto the basis set) so subtraction, just like that of "standard"
|
Chris@87
|
378 polynomials, is simply "component-wise."
|
Chris@87
|
379
|
Chris@87
|
380 Examples
|
Chris@87
|
381 --------
|
Chris@87
|
382 >>> from numpy.polynomial.laguerre import lagsub
|
Chris@87
|
383 >>> lagsub([1, 2, 3, 4], [1, 2, 3])
|
Chris@87
|
384 array([ 0., 0., 0., 4.])
|
Chris@87
|
385
|
Chris@87
|
386 """
|
Chris@87
|
387 # c1, c2 are trimmed copies
|
Chris@87
|
388 [c1, c2] = pu.as_series([c1, c2])
|
Chris@87
|
389 if len(c1) > len(c2):
|
Chris@87
|
390 c1[:c2.size] -= c2
|
Chris@87
|
391 ret = c1
|
Chris@87
|
392 else:
|
Chris@87
|
393 c2 = -c2
|
Chris@87
|
394 c2[:c1.size] += c1
|
Chris@87
|
395 ret = c2
|
Chris@87
|
396 return pu.trimseq(ret)
|
Chris@87
|
397
|
Chris@87
|
398
|
Chris@87
|
399 def lagmulx(c):
|
Chris@87
|
400 """Multiply a Laguerre series by x.
|
Chris@87
|
401
|
Chris@87
|
402 Multiply the Laguerre series `c` by x, where x is the independent
|
Chris@87
|
403 variable.
|
Chris@87
|
404
|
Chris@87
|
405
|
Chris@87
|
406 Parameters
|
Chris@87
|
407 ----------
|
Chris@87
|
408 c : array_like
|
Chris@87
|
409 1-D array of Laguerre series coefficients ordered from low to
|
Chris@87
|
410 high.
|
Chris@87
|
411
|
Chris@87
|
412 Returns
|
Chris@87
|
413 -------
|
Chris@87
|
414 out : ndarray
|
Chris@87
|
415 Array representing the result of the multiplication.
|
Chris@87
|
416
|
Chris@87
|
417 Notes
|
Chris@87
|
418 -----
|
Chris@87
|
419 The multiplication uses the recursion relationship for Laguerre
|
Chris@87
|
420 polynomials in the form
|
Chris@87
|
421
|
Chris@87
|
422 .. math::
|
Chris@87
|
423
|
Chris@87
|
424 xP_i(x) = (-(i + 1)*P_{i + 1}(x) + (2i + 1)P_{i}(x) - iP_{i - 1}(x))
|
Chris@87
|
425
|
Chris@87
|
426 Examples
|
Chris@87
|
427 --------
|
Chris@87
|
428 >>> from numpy.polynomial.laguerre import lagmulx
|
Chris@87
|
429 >>> lagmulx([1, 2, 3])
|
Chris@87
|
430 array([ -1., -1., 11., -9.])
|
Chris@87
|
431
|
Chris@87
|
432 """
|
Chris@87
|
433 # c is a trimmed copy
|
Chris@87
|
434 [c] = pu.as_series([c])
|
Chris@87
|
435 # The zero series needs special treatment
|
Chris@87
|
436 if len(c) == 1 and c[0] == 0:
|
Chris@87
|
437 return c
|
Chris@87
|
438
|
Chris@87
|
439 prd = np.empty(len(c) + 1, dtype=c.dtype)
|
Chris@87
|
440 prd[0] = c[0]
|
Chris@87
|
441 prd[1] = -c[0]
|
Chris@87
|
442 for i in range(1, len(c)):
|
Chris@87
|
443 prd[i + 1] = -c[i]*(i + 1)
|
Chris@87
|
444 prd[i] += c[i]*(2*i + 1)
|
Chris@87
|
445 prd[i - 1] -= c[i]*i
|
Chris@87
|
446 return prd
|
Chris@87
|
447
|
Chris@87
|
448
|
Chris@87
|
449 def lagmul(c1, c2):
|
Chris@87
|
450 """
|
Chris@87
|
451 Multiply one Laguerre series by another.
|
Chris@87
|
452
|
Chris@87
|
453 Returns the product of two Laguerre series `c1` * `c2`. The arguments
|
Chris@87
|
454 are sequences of coefficients, from lowest order "term" to highest,
|
Chris@87
|
455 e.g., [1,2,3] represents the series ``P_0 + 2*P_1 + 3*P_2``.
|
Chris@87
|
456
|
Chris@87
|
457 Parameters
|
Chris@87
|
458 ----------
|
Chris@87
|
459 c1, c2 : array_like
|
Chris@87
|
460 1-D arrays of Laguerre series coefficients ordered from low to
|
Chris@87
|
461 high.
|
Chris@87
|
462
|
Chris@87
|
463 Returns
|
Chris@87
|
464 -------
|
Chris@87
|
465 out : ndarray
|
Chris@87
|
466 Of Laguerre series coefficients representing their product.
|
Chris@87
|
467
|
Chris@87
|
468 See Also
|
Chris@87
|
469 --------
|
Chris@87
|
470 lagadd, lagsub, lagdiv, lagpow
|
Chris@87
|
471
|
Chris@87
|
472 Notes
|
Chris@87
|
473 -----
|
Chris@87
|
474 In general, the (polynomial) product of two C-series results in terms
|
Chris@87
|
475 that are not in the Laguerre polynomial basis set. Thus, to express
|
Chris@87
|
476 the product as a Laguerre series, it is necessary to "reproject" the
|
Chris@87
|
477 product onto said basis set, which may produce "unintuitive" (but
|
Chris@87
|
478 correct) results; see Examples section below.
|
Chris@87
|
479
|
Chris@87
|
480 Examples
|
Chris@87
|
481 --------
|
Chris@87
|
482 >>> from numpy.polynomial.laguerre import lagmul
|
Chris@87
|
483 >>> lagmul([1, 2, 3], [0, 1, 2])
|
Chris@87
|
484 array([ 8., -13., 38., -51., 36.])
|
Chris@87
|
485
|
Chris@87
|
486 """
|
Chris@87
|
487 # s1, s2 are trimmed copies
|
Chris@87
|
488 [c1, c2] = pu.as_series([c1, c2])
|
Chris@87
|
489
|
Chris@87
|
490 if len(c1) > len(c2):
|
Chris@87
|
491 c = c2
|
Chris@87
|
492 xs = c1
|
Chris@87
|
493 else:
|
Chris@87
|
494 c = c1
|
Chris@87
|
495 xs = c2
|
Chris@87
|
496
|
Chris@87
|
497 if len(c) == 1:
|
Chris@87
|
498 c0 = c[0]*xs
|
Chris@87
|
499 c1 = 0
|
Chris@87
|
500 elif len(c) == 2:
|
Chris@87
|
501 c0 = c[0]*xs
|
Chris@87
|
502 c1 = c[1]*xs
|
Chris@87
|
503 else:
|
Chris@87
|
504 nd = len(c)
|
Chris@87
|
505 c0 = c[-2]*xs
|
Chris@87
|
506 c1 = c[-1]*xs
|
Chris@87
|
507 for i in range(3, len(c) + 1):
|
Chris@87
|
508 tmp = c0
|
Chris@87
|
509 nd = nd - 1
|
Chris@87
|
510 c0 = lagsub(c[-i]*xs, (c1*(nd - 1))/nd)
|
Chris@87
|
511 c1 = lagadd(tmp, lagsub((2*nd - 1)*c1, lagmulx(c1))/nd)
|
Chris@87
|
512 return lagadd(c0, lagsub(c1, lagmulx(c1)))
|
Chris@87
|
513
|
Chris@87
|
514
|
Chris@87
|
515 def lagdiv(c1, c2):
|
Chris@87
|
516 """
|
Chris@87
|
517 Divide one Laguerre series by another.
|
Chris@87
|
518
|
Chris@87
|
519 Returns the quotient-with-remainder of two Laguerre series
|
Chris@87
|
520 `c1` / `c2`. The arguments are sequences of coefficients from lowest
|
Chris@87
|
521 order "term" to highest, e.g., [1,2,3] represents the series
|
Chris@87
|
522 ``P_0 + 2*P_1 + 3*P_2``.
|
Chris@87
|
523
|
Chris@87
|
524 Parameters
|
Chris@87
|
525 ----------
|
Chris@87
|
526 c1, c2 : array_like
|
Chris@87
|
527 1-D arrays of Laguerre series coefficients ordered from low to
|
Chris@87
|
528 high.
|
Chris@87
|
529
|
Chris@87
|
530 Returns
|
Chris@87
|
531 -------
|
Chris@87
|
532 [quo, rem] : ndarrays
|
Chris@87
|
533 Of Laguerre series coefficients representing the quotient and
|
Chris@87
|
534 remainder.
|
Chris@87
|
535
|
Chris@87
|
536 See Also
|
Chris@87
|
537 --------
|
Chris@87
|
538 lagadd, lagsub, lagmul, lagpow
|
Chris@87
|
539
|
Chris@87
|
540 Notes
|
Chris@87
|
541 -----
|
Chris@87
|
542 In general, the (polynomial) division of one Laguerre series by another
|
Chris@87
|
543 results in quotient and remainder terms that are not in the Laguerre
|
Chris@87
|
544 polynomial basis set. Thus, to express these results as a Laguerre
|
Chris@87
|
545 series, it is necessary to "reproject" the results onto the Laguerre
|
Chris@87
|
546 basis set, which may produce "unintuitive" (but correct) results; see
|
Chris@87
|
547 Examples section below.
|
Chris@87
|
548
|
Chris@87
|
549 Examples
|
Chris@87
|
550 --------
|
Chris@87
|
551 >>> from numpy.polynomial.laguerre import lagdiv
|
Chris@87
|
552 >>> lagdiv([ 8., -13., 38., -51., 36.], [0, 1, 2])
|
Chris@87
|
553 (array([ 1., 2., 3.]), array([ 0.]))
|
Chris@87
|
554 >>> lagdiv([ 9., -12., 38., -51., 36.], [0, 1, 2])
|
Chris@87
|
555 (array([ 1., 2., 3.]), array([ 1., 1.]))
|
Chris@87
|
556
|
Chris@87
|
557 """
|
Chris@87
|
558 # c1, c2 are trimmed copies
|
Chris@87
|
559 [c1, c2] = pu.as_series([c1, c2])
|
Chris@87
|
560 if c2[-1] == 0:
|
Chris@87
|
561 raise ZeroDivisionError()
|
Chris@87
|
562
|
Chris@87
|
563 lc1 = len(c1)
|
Chris@87
|
564 lc2 = len(c2)
|
Chris@87
|
565 if lc1 < lc2:
|
Chris@87
|
566 return c1[:1]*0, c1
|
Chris@87
|
567 elif lc2 == 1:
|
Chris@87
|
568 return c1/c2[-1], c1[:1]*0
|
Chris@87
|
569 else:
|
Chris@87
|
570 quo = np.empty(lc1 - lc2 + 1, dtype=c1.dtype)
|
Chris@87
|
571 rem = c1
|
Chris@87
|
572 for i in range(lc1 - lc2, - 1, -1):
|
Chris@87
|
573 p = lagmul([0]*i + [1], c2)
|
Chris@87
|
574 q = rem[-1]/p[-1]
|
Chris@87
|
575 rem = rem[:-1] - q*p[:-1]
|
Chris@87
|
576 quo[i] = q
|
Chris@87
|
577 return quo, pu.trimseq(rem)
|
Chris@87
|
578
|
Chris@87
|
579
|
Chris@87
|
580 def lagpow(c, pow, maxpower=16):
|
Chris@87
|
581 """Raise a Laguerre series to a power.
|
Chris@87
|
582
|
Chris@87
|
583 Returns the Laguerre series `c` raised to the power `pow`. The
|
Chris@87
|
584 argument `c` is a sequence of coefficients ordered from low to high.
|
Chris@87
|
585 i.e., [1,2,3] is the series ``P_0 + 2*P_1 + 3*P_2.``
|
Chris@87
|
586
|
Chris@87
|
587 Parameters
|
Chris@87
|
588 ----------
|
Chris@87
|
589 c : array_like
|
Chris@87
|
590 1-D array of Laguerre series coefficients ordered from low to
|
Chris@87
|
591 high.
|
Chris@87
|
592 pow : integer
|
Chris@87
|
593 Power to which the series will be raised
|
Chris@87
|
594 maxpower : integer, optional
|
Chris@87
|
595 Maximum power allowed. This is mainly to limit growth of the series
|
Chris@87
|
596 to unmanageable size. Default is 16
|
Chris@87
|
597
|
Chris@87
|
598 Returns
|
Chris@87
|
599 -------
|
Chris@87
|
600 coef : ndarray
|
Chris@87
|
601 Laguerre series of power.
|
Chris@87
|
602
|
Chris@87
|
603 See Also
|
Chris@87
|
604 --------
|
Chris@87
|
605 lagadd, lagsub, lagmul, lagdiv
|
Chris@87
|
606
|
Chris@87
|
607 Examples
|
Chris@87
|
608 --------
|
Chris@87
|
609 >>> from numpy.polynomial.laguerre import lagpow
|
Chris@87
|
610 >>> lagpow([1, 2, 3], 2)
|
Chris@87
|
611 array([ 14., -16., 56., -72., 54.])
|
Chris@87
|
612
|
Chris@87
|
613 """
|
Chris@87
|
614 # c is a trimmed copy
|
Chris@87
|
615 [c] = pu.as_series([c])
|
Chris@87
|
616 power = int(pow)
|
Chris@87
|
617 if power != pow or power < 0:
|
Chris@87
|
618 raise ValueError("Power must be a non-negative integer.")
|
Chris@87
|
619 elif maxpower is not None and power > maxpower:
|
Chris@87
|
620 raise ValueError("Power is too large")
|
Chris@87
|
621 elif power == 0:
|
Chris@87
|
622 return np.array([1], dtype=c.dtype)
|
Chris@87
|
623 elif power == 1:
|
Chris@87
|
624 return c
|
Chris@87
|
625 else:
|
Chris@87
|
626 # This can be made more efficient by using powers of two
|
Chris@87
|
627 # in the usual way.
|
Chris@87
|
628 prd = c
|
Chris@87
|
629 for i in range(2, power + 1):
|
Chris@87
|
630 prd = lagmul(prd, c)
|
Chris@87
|
631 return prd
|
Chris@87
|
632
|
Chris@87
|
633
|
Chris@87
|
634 def lagder(c, m=1, scl=1, axis=0):
|
Chris@87
|
635 """
|
Chris@87
|
636 Differentiate a Laguerre series.
|
Chris@87
|
637
|
Chris@87
|
638 Returns the Laguerre series coefficients `c` differentiated `m` times
|
Chris@87
|
639 along `axis`. At each iteration the result is multiplied by `scl` (the
|
Chris@87
|
640 scaling factor is for use in a linear change of variable). The argument
|
Chris@87
|
641 `c` is an array of coefficients from low to high degree along each
|
Chris@87
|
642 axis, e.g., [1,2,3] represents the series ``1*L_0 + 2*L_1 + 3*L_2``
|
Chris@87
|
643 while [[1,2],[1,2]] represents ``1*L_0(x)*L_0(y) + 1*L_1(x)*L_0(y) +
|
Chris@87
|
644 2*L_0(x)*L_1(y) + 2*L_1(x)*L_1(y)`` if axis=0 is ``x`` and axis=1 is
|
Chris@87
|
645 ``y``.
|
Chris@87
|
646
|
Chris@87
|
647 Parameters
|
Chris@87
|
648 ----------
|
Chris@87
|
649 c : array_like
|
Chris@87
|
650 Array of Laguerre series coefficients. If `c` is multidimensional
|
Chris@87
|
651 the different axis correspond to different variables with the
|
Chris@87
|
652 degree in each axis given by the corresponding index.
|
Chris@87
|
653 m : int, optional
|
Chris@87
|
654 Number of derivatives taken, must be non-negative. (Default: 1)
|
Chris@87
|
655 scl : scalar, optional
|
Chris@87
|
656 Each differentiation is multiplied by `scl`. The end result is
|
Chris@87
|
657 multiplication by ``scl**m``. This is for use in a linear change of
|
Chris@87
|
658 variable. (Default: 1)
|
Chris@87
|
659 axis : int, optional
|
Chris@87
|
660 Axis over which the derivative is taken. (Default: 0).
|
Chris@87
|
661
|
Chris@87
|
662 .. versionadded:: 1.7.0
|
Chris@87
|
663
|
Chris@87
|
664 Returns
|
Chris@87
|
665 -------
|
Chris@87
|
666 der : ndarray
|
Chris@87
|
667 Laguerre series of the derivative.
|
Chris@87
|
668
|
Chris@87
|
669 See Also
|
Chris@87
|
670 --------
|
Chris@87
|
671 lagint
|
Chris@87
|
672
|
Chris@87
|
673 Notes
|
Chris@87
|
674 -----
|
Chris@87
|
675 In general, the result of differentiating a Laguerre series does not
|
Chris@87
|
676 resemble the same operation on a power series. Thus the result of this
|
Chris@87
|
677 function may be "unintuitive," albeit correct; see Examples section
|
Chris@87
|
678 below.
|
Chris@87
|
679
|
Chris@87
|
680 Examples
|
Chris@87
|
681 --------
|
Chris@87
|
682 >>> from numpy.polynomial.laguerre import lagder
|
Chris@87
|
683 >>> lagder([ 1., 1., 1., -3.])
|
Chris@87
|
684 array([ 1., 2., 3.])
|
Chris@87
|
685 >>> lagder([ 1., 0., 0., -4., 3.], m=2)
|
Chris@87
|
686 array([ 1., 2., 3.])
|
Chris@87
|
687
|
Chris@87
|
688 """
|
Chris@87
|
689 c = np.array(c, ndmin=1, copy=1)
|
Chris@87
|
690 if c.dtype.char in '?bBhHiIlLqQpP':
|
Chris@87
|
691 c = c.astype(np.double)
|
Chris@87
|
692 cnt, iaxis = [int(t) for t in [m, axis]]
|
Chris@87
|
693
|
Chris@87
|
694 if cnt != m:
|
Chris@87
|
695 raise ValueError("The order of derivation must be integer")
|
Chris@87
|
696 if cnt < 0:
|
Chris@87
|
697 raise ValueError("The order of derivation must be non-negative")
|
Chris@87
|
698 if iaxis != axis:
|
Chris@87
|
699 raise ValueError("The axis must be integer")
|
Chris@87
|
700 if not -c.ndim <= iaxis < c.ndim:
|
Chris@87
|
701 raise ValueError("The axis is out of range")
|
Chris@87
|
702 if iaxis < 0:
|
Chris@87
|
703 iaxis += c.ndim
|
Chris@87
|
704
|
Chris@87
|
705 if cnt == 0:
|
Chris@87
|
706 return c
|
Chris@87
|
707
|
Chris@87
|
708 c = np.rollaxis(c, iaxis)
|
Chris@87
|
709 n = len(c)
|
Chris@87
|
710 if cnt >= n:
|
Chris@87
|
711 c = c[:1]*0
|
Chris@87
|
712 else:
|
Chris@87
|
713 for i in range(cnt):
|
Chris@87
|
714 n = n - 1
|
Chris@87
|
715 c *= scl
|
Chris@87
|
716 der = np.empty((n,) + c.shape[1:], dtype=c.dtype)
|
Chris@87
|
717 for j in range(n, 1, -1):
|
Chris@87
|
718 der[j - 1] = -c[j]
|
Chris@87
|
719 c[j - 1] += c[j]
|
Chris@87
|
720 der[0] = -c[1]
|
Chris@87
|
721 c = der
|
Chris@87
|
722 c = np.rollaxis(c, 0, iaxis + 1)
|
Chris@87
|
723 return c
|
Chris@87
|
724
|
Chris@87
|
725
|
Chris@87
|
726 def lagint(c, m=1, k=[], lbnd=0, scl=1, axis=0):
|
Chris@87
|
727 """
|
Chris@87
|
728 Integrate a Laguerre series.
|
Chris@87
|
729
|
Chris@87
|
730 Returns the Laguerre series coefficients `c` integrated `m` times from
|
Chris@87
|
731 `lbnd` along `axis`. At each iteration the resulting series is
|
Chris@87
|
732 **multiplied** by `scl` and an integration constant, `k`, is added.
|
Chris@87
|
733 The scaling factor is for use in a linear change of variable. ("Buyer
|
Chris@87
|
734 beware": note that, depending on what one is doing, one may want `scl`
|
Chris@87
|
735 to be the reciprocal of what one might expect; for more information,
|
Chris@87
|
736 see the Notes section below.) The argument `c` is an array of
|
Chris@87
|
737 coefficients from low to high degree along each axis, e.g., [1,2,3]
|
Chris@87
|
738 represents the series ``L_0 + 2*L_1 + 3*L_2`` while [[1,2],[1,2]]
|
Chris@87
|
739 represents ``1*L_0(x)*L_0(y) + 1*L_1(x)*L_0(y) + 2*L_0(x)*L_1(y) +
|
Chris@87
|
740 2*L_1(x)*L_1(y)`` if axis=0 is ``x`` and axis=1 is ``y``.
|
Chris@87
|
741
|
Chris@87
|
742
|
Chris@87
|
743 Parameters
|
Chris@87
|
744 ----------
|
Chris@87
|
745 c : array_like
|
Chris@87
|
746 Array of Laguerre series coefficients. If `c` is multidimensional
|
Chris@87
|
747 the different axis correspond to different variables with the
|
Chris@87
|
748 degree in each axis given by the corresponding index.
|
Chris@87
|
749 m : int, optional
|
Chris@87
|
750 Order of integration, must be positive. (Default: 1)
|
Chris@87
|
751 k : {[], list, scalar}, optional
|
Chris@87
|
752 Integration constant(s). The value of the first integral at
|
Chris@87
|
753 ``lbnd`` is the first value in the list, the value of the second
|
Chris@87
|
754 integral at ``lbnd`` is the second value, etc. If ``k == []`` (the
|
Chris@87
|
755 default), all constants are set to zero. If ``m == 1``, a single
|
Chris@87
|
756 scalar can be given instead of a list.
|
Chris@87
|
757 lbnd : scalar, optional
|
Chris@87
|
758 The lower bound of the integral. (Default: 0)
|
Chris@87
|
759 scl : scalar, optional
|
Chris@87
|
760 Following each integration the result is *multiplied* by `scl`
|
Chris@87
|
761 before the integration constant is added. (Default: 1)
|
Chris@87
|
762 axis : int, optional
|
Chris@87
|
763 Axis over which the integral is taken. (Default: 0).
|
Chris@87
|
764
|
Chris@87
|
765 .. versionadded:: 1.7.0
|
Chris@87
|
766
|
Chris@87
|
767 Returns
|
Chris@87
|
768 -------
|
Chris@87
|
769 S : ndarray
|
Chris@87
|
770 Laguerre series coefficients of the integral.
|
Chris@87
|
771
|
Chris@87
|
772 Raises
|
Chris@87
|
773 ------
|
Chris@87
|
774 ValueError
|
Chris@87
|
775 If ``m < 0``, ``len(k) > m``, ``np.isscalar(lbnd) == False``, or
|
Chris@87
|
776 ``np.isscalar(scl) == False``.
|
Chris@87
|
777
|
Chris@87
|
778 See Also
|
Chris@87
|
779 --------
|
Chris@87
|
780 lagder
|
Chris@87
|
781
|
Chris@87
|
782 Notes
|
Chris@87
|
783 -----
|
Chris@87
|
784 Note that the result of each integration is *multiplied* by `scl`.
|
Chris@87
|
785 Why is this important to note? Say one is making a linear change of
|
Chris@87
|
786 variable :math:`u = ax + b` in an integral relative to `x`. Then
|
Chris@87
|
787 .. math::`dx = du/a`, so one will need to set `scl` equal to
|
Chris@87
|
788 :math:`1/a` - perhaps not what one would have first thought.
|
Chris@87
|
789
|
Chris@87
|
790 Also note that, in general, the result of integrating a C-series needs
|
Chris@87
|
791 to be "reprojected" onto the C-series basis set. Thus, typically,
|
Chris@87
|
792 the result of this function is "unintuitive," albeit correct; see
|
Chris@87
|
793 Examples section below.
|
Chris@87
|
794
|
Chris@87
|
795 Examples
|
Chris@87
|
796 --------
|
Chris@87
|
797 >>> from numpy.polynomial.laguerre import lagint
|
Chris@87
|
798 >>> lagint([1,2,3])
|
Chris@87
|
799 array([ 1., 1., 1., -3.])
|
Chris@87
|
800 >>> lagint([1,2,3], m=2)
|
Chris@87
|
801 array([ 1., 0., 0., -4., 3.])
|
Chris@87
|
802 >>> lagint([1,2,3], k=1)
|
Chris@87
|
803 array([ 2., 1., 1., -3.])
|
Chris@87
|
804 >>> lagint([1,2,3], lbnd=-1)
|
Chris@87
|
805 array([ 11.5, 1. , 1. , -3. ])
|
Chris@87
|
806 >>> lagint([1,2], m=2, k=[1,2], lbnd=-1)
|
Chris@87
|
807 array([ 11.16666667, -5. , -3. , 2. ])
|
Chris@87
|
808
|
Chris@87
|
809 """
|
Chris@87
|
810 c = np.array(c, ndmin=1, copy=1)
|
Chris@87
|
811 if c.dtype.char in '?bBhHiIlLqQpP':
|
Chris@87
|
812 c = c.astype(np.double)
|
Chris@87
|
813 if not np.iterable(k):
|
Chris@87
|
814 k = [k]
|
Chris@87
|
815 cnt, iaxis = [int(t) for t in [m, axis]]
|
Chris@87
|
816
|
Chris@87
|
817 if cnt != m:
|
Chris@87
|
818 raise ValueError("The order of integration must be integer")
|
Chris@87
|
819 if cnt < 0:
|
Chris@87
|
820 raise ValueError("The order of integration must be non-negative")
|
Chris@87
|
821 if len(k) > cnt:
|
Chris@87
|
822 raise ValueError("Too many integration constants")
|
Chris@87
|
823 if iaxis != axis:
|
Chris@87
|
824 raise ValueError("The axis must be integer")
|
Chris@87
|
825 if not -c.ndim <= iaxis < c.ndim:
|
Chris@87
|
826 raise ValueError("The axis is out of range")
|
Chris@87
|
827 if iaxis < 0:
|
Chris@87
|
828 iaxis += c.ndim
|
Chris@87
|
829
|
Chris@87
|
830 if cnt == 0:
|
Chris@87
|
831 return c
|
Chris@87
|
832
|
Chris@87
|
833 c = np.rollaxis(c, iaxis)
|
Chris@87
|
834 k = list(k) + [0]*(cnt - len(k))
|
Chris@87
|
835 for i in range(cnt):
|
Chris@87
|
836 n = len(c)
|
Chris@87
|
837 c *= scl
|
Chris@87
|
838 if n == 1 and np.all(c[0] == 0):
|
Chris@87
|
839 c[0] += k[i]
|
Chris@87
|
840 else:
|
Chris@87
|
841 tmp = np.empty((n + 1,) + c.shape[1:], dtype=c.dtype)
|
Chris@87
|
842 tmp[0] = c[0]
|
Chris@87
|
843 tmp[1] = -c[0]
|
Chris@87
|
844 for j in range(1, n):
|
Chris@87
|
845 tmp[j] += c[j]
|
Chris@87
|
846 tmp[j + 1] = -c[j]
|
Chris@87
|
847 tmp[0] += k[i] - lagval(lbnd, tmp)
|
Chris@87
|
848 c = tmp
|
Chris@87
|
849 c = np.rollaxis(c, 0, iaxis + 1)
|
Chris@87
|
850 return c
|
Chris@87
|
851
|
Chris@87
|
852
|
Chris@87
|
853 def lagval(x, c, tensor=True):
|
Chris@87
|
854 """
|
Chris@87
|
855 Evaluate a Laguerre series at points x.
|
Chris@87
|
856
|
Chris@87
|
857 If `c` is of length `n + 1`, this function returns the value:
|
Chris@87
|
858
|
Chris@87
|
859 .. math:: p(x) = c_0 * L_0(x) + c_1 * L_1(x) + ... + c_n * L_n(x)
|
Chris@87
|
860
|
Chris@87
|
861 The parameter `x` is converted to an array only if it is a tuple or a
|
Chris@87
|
862 list, otherwise it is treated as a scalar. In either case, either `x`
|
Chris@87
|
863 or its elements must support multiplication and addition both with
|
Chris@87
|
864 themselves and with the elements of `c`.
|
Chris@87
|
865
|
Chris@87
|
866 If `c` is a 1-D array, then `p(x)` will have the same shape as `x`. If
|
Chris@87
|
867 `c` is multidimensional, then the shape of the result depends on the
|
Chris@87
|
868 value of `tensor`. If `tensor` is true the shape will be c.shape[1:] +
|
Chris@87
|
869 x.shape. If `tensor` is false the shape will be c.shape[1:]. Note that
|
Chris@87
|
870 scalars have shape (,).
|
Chris@87
|
871
|
Chris@87
|
872 Trailing zeros in the coefficients will be used in the evaluation, so
|
Chris@87
|
873 they should be avoided if efficiency is a concern.
|
Chris@87
|
874
|
Chris@87
|
875 Parameters
|
Chris@87
|
876 ----------
|
Chris@87
|
877 x : array_like, compatible object
|
Chris@87
|
878 If `x` is a list or tuple, it is converted to an ndarray, otherwise
|
Chris@87
|
879 it is left unchanged and treated as a scalar. In either case, `x`
|
Chris@87
|
880 or its elements must support addition and multiplication with
|
Chris@87
|
881 with themselves and with the elements of `c`.
|
Chris@87
|
882 c : array_like
|
Chris@87
|
883 Array of coefficients ordered so that the coefficients for terms of
|
Chris@87
|
884 degree n are contained in c[n]. If `c` is multidimensional the
|
Chris@87
|
885 remaining indices enumerate multiple polynomials. In the two
|
Chris@87
|
886 dimensional case the coefficients may be thought of as stored in
|
Chris@87
|
887 the columns of `c`.
|
Chris@87
|
888 tensor : boolean, optional
|
Chris@87
|
889 If True, the shape of the coefficient array is extended with ones
|
Chris@87
|
890 on the right, one for each dimension of `x`. Scalars have dimension 0
|
Chris@87
|
891 for this action. The result is that every column of coefficients in
|
Chris@87
|
892 `c` is evaluated for every element of `x`. If False, `x` is broadcast
|
Chris@87
|
893 over the columns of `c` for the evaluation. This keyword is useful
|
Chris@87
|
894 when `c` is multidimensional. The default value is True.
|
Chris@87
|
895
|
Chris@87
|
896 .. versionadded:: 1.7.0
|
Chris@87
|
897
|
Chris@87
|
898 Returns
|
Chris@87
|
899 -------
|
Chris@87
|
900 values : ndarray, algebra_like
|
Chris@87
|
901 The shape of the return value is described above.
|
Chris@87
|
902
|
Chris@87
|
903 See Also
|
Chris@87
|
904 --------
|
Chris@87
|
905 lagval2d, laggrid2d, lagval3d, laggrid3d
|
Chris@87
|
906
|
Chris@87
|
907 Notes
|
Chris@87
|
908 -----
|
Chris@87
|
909 The evaluation uses Clenshaw recursion, aka synthetic division.
|
Chris@87
|
910
|
Chris@87
|
911 Examples
|
Chris@87
|
912 --------
|
Chris@87
|
913 >>> from numpy.polynomial.laguerre import lagval
|
Chris@87
|
914 >>> coef = [1,2,3]
|
Chris@87
|
915 >>> lagval(1, coef)
|
Chris@87
|
916 -0.5
|
Chris@87
|
917 >>> lagval([[1,2],[3,4]], coef)
|
Chris@87
|
918 array([[-0.5, -4. ],
|
Chris@87
|
919 [-4.5, -2. ]])
|
Chris@87
|
920
|
Chris@87
|
921 """
|
Chris@87
|
922 c = np.array(c, ndmin=1, copy=0)
|
Chris@87
|
923 if c.dtype.char in '?bBhHiIlLqQpP':
|
Chris@87
|
924 c = c.astype(np.double)
|
Chris@87
|
925 if isinstance(x, (tuple, list)):
|
Chris@87
|
926 x = np.asarray(x)
|
Chris@87
|
927 if isinstance(x, np.ndarray) and tensor:
|
Chris@87
|
928 c = c.reshape(c.shape + (1,)*x.ndim)
|
Chris@87
|
929
|
Chris@87
|
930 if len(c) == 1:
|
Chris@87
|
931 c0 = c[0]
|
Chris@87
|
932 c1 = 0
|
Chris@87
|
933 elif len(c) == 2:
|
Chris@87
|
934 c0 = c[0]
|
Chris@87
|
935 c1 = c[1]
|
Chris@87
|
936 else:
|
Chris@87
|
937 nd = len(c)
|
Chris@87
|
938 c0 = c[-2]
|
Chris@87
|
939 c1 = c[-1]
|
Chris@87
|
940 for i in range(3, len(c) + 1):
|
Chris@87
|
941 tmp = c0
|
Chris@87
|
942 nd = nd - 1
|
Chris@87
|
943 c0 = c[-i] - (c1*(nd - 1))/nd
|
Chris@87
|
944 c1 = tmp + (c1*((2*nd - 1) - x))/nd
|
Chris@87
|
945 return c0 + c1*(1 - x)
|
Chris@87
|
946
|
Chris@87
|
947
|
Chris@87
|
948 def lagval2d(x, y, c):
|
Chris@87
|
949 """
|
Chris@87
|
950 Evaluate a 2-D Laguerre series at points (x, y).
|
Chris@87
|
951
|
Chris@87
|
952 This function returns the values:
|
Chris@87
|
953
|
Chris@87
|
954 .. math:: p(x,y) = \\sum_{i,j} c_{i,j} * L_i(x) * L_j(y)
|
Chris@87
|
955
|
Chris@87
|
956 The parameters `x` and `y` are converted to arrays only if they are
|
Chris@87
|
957 tuples or a lists, otherwise they are treated as a scalars and they
|
Chris@87
|
958 must have the same shape after conversion. In either case, either `x`
|
Chris@87
|
959 and `y` or their elements must support multiplication and addition both
|
Chris@87
|
960 with themselves and with the elements of `c`.
|
Chris@87
|
961
|
Chris@87
|
962 If `c` is a 1-D array a one is implicitly appended to its shape to make
|
Chris@87
|
963 it 2-D. The shape of the result will be c.shape[2:] + x.shape.
|
Chris@87
|
964
|
Chris@87
|
965 Parameters
|
Chris@87
|
966 ----------
|
Chris@87
|
967 x, y : array_like, compatible objects
|
Chris@87
|
968 The two dimensional series is evaluated at the points `(x, y)`,
|
Chris@87
|
969 where `x` and `y` must have the same shape. If `x` or `y` is a list
|
Chris@87
|
970 or tuple, it is first converted to an ndarray, otherwise it is left
|
Chris@87
|
971 unchanged and if it isn't an ndarray it is treated as a scalar.
|
Chris@87
|
972 c : array_like
|
Chris@87
|
973 Array of coefficients ordered so that the coefficient of the term
|
Chris@87
|
974 of multi-degree i,j is contained in ``c[i,j]``. If `c` has
|
Chris@87
|
975 dimension greater than two the remaining indices enumerate multiple
|
Chris@87
|
976 sets of coefficients.
|
Chris@87
|
977
|
Chris@87
|
978 Returns
|
Chris@87
|
979 -------
|
Chris@87
|
980 values : ndarray, compatible object
|
Chris@87
|
981 The values of the two dimensional polynomial at points formed with
|
Chris@87
|
982 pairs of corresponding values from `x` and `y`.
|
Chris@87
|
983
|
Chris@87
|
984 See Also
|
Chris@87
|
985 --------
|
Chris@87
|
986 lagval, laggrid2d, lagval3d, laggrid3d
|
Chris@87
|
987
|
Chris@87
|
988 Notes
|
Chris@87
|
989 -----
|
Chris@87
|
990
|
Chris@87
|
991 .. versionadded::1.7.0
|
Chris@87
|
992
|
Chris@87
|
993 """
|
Chris@87
|
994 try:
|
Chris@87
|
995 x, y = np.array((x, y), copy=0)
|
Chris@87
|
996 except:
|
Chris@87
|
997 raise ValueError('x, y are incompatible')
|
Chris@87
|
998
|
Chris@87
|
999 c = lagval(x, c)
|
Chris@87
|
1000 c = lagval(y, c, tensor=False)
|
Chris@87
|
1001 return c
|
Chris@87
|
1002
|
Chris@87
|
1003
|
Chris@87
|
1004 def laggrid2d(x, y, c):
|
Chris@87
|
1005 """
|
Chris@87
|
1006 Evaluate a 2-D Laguerre series on the Cartesian product of x and y.
|
Chris@87
|
1007
|
Chris@87
|
1008 This function returns the values:
|
Chris@87
|
1009
|
Chris@87
|
1010 .. math:: p(a,b) = \sum_{i,j} c_{i,j} * L_i(a) * L_j(b)
|
Chris@87
|
1011
|
Chris@87
|
1012 where the points `(a, b)` consist of all pairs formed by taking
|
Chris@87
|
1013 `a` from `x` and `b` from `y`. The resulting points form a grid with
|
Chris@87
|
1014 `x` in the first dimension and `y` in the second.
|
Chris@87
|
1015
|
Chris@87
|
1016 The parameters `x` and `y` are converted to arrays only if they are
|
Chris@87
|
1017 tuples or a lists, otherwise they are treated as a scalars. In either
|
Chris@87
|
1018 case, either `x` and `y` or their elements must support multiplication
|
Chris@87
|
1019 and addition both with themselves and with the elements of `c`.
|
Chris@87
|
1020
|
Chris@87
|
1021 If `c` has fewer than two dimensions, ones are implicitly appended to
|
Chris@87
|
1022 its shape to make it 2-D. The shape of the result will be c.shape[2:] +
|
Chris@87
|
1023 x.shape + y.shape.
|
Chris@87
|
1024
|
Chris@87
|
1025 Parameters
|
Chris@87
|
1026 ----------
|
Chris@87
|
1027 x, y : array_like, compatible objects
|
Chris@87
|
1028 The two dimensional series is evaluated at the points in the
|
Chris@87
|
1029 Cartesian product of `x` and `y`. If `x` or `y` is a list or
|
Chris@87
|
1030 tuple, it is first converted to an ndarray, otherwise it is left
|
Chris@87
|
1031 unchanged and, if it isn't an ndarray, it is treated as a scalar.
|
Chris@87
|
1032 c : array_like
|
Chris@87
|
1033 Array of coefficients ordered so that the coefficient of the term of
|
Chris@87
|
1034 multi-degree i,j is contained in `c[i,j]`. If `c` has dimension
|
Chris@87
|
1035 greater than two the remaining indices enumerate multiple sets of
|
Chris@87
|
1036 coefficients.
|
Chris@87
|
1037
|
Chris@87
|
1038 Returns
|
Chris@87
|
1039 -------
|
Chris@87
|
1040 values : ndarray, compatible object
|
Chris@87
|
1041 The values of the two dimensional Chebyshev series at points in the
|
Chris@87
|
1042 Cartesian product of `x` and `y`.
|
Chris@87
|
1043
|
Chris@87
|
1044 See Also
|
Chris@87
|
1045 --------
|
Chris@87
|
1046 lagval, lagval2d, lagval3d, laggrid3d
|
Chris@87
|
1047
|
Chris@87
|
1048 Notes
|
Chris@87
|
1049 -----
|
Chris@87
|
1050
|
Chris@87
|
1051 .. versionadded::1.7.0
|
Chris@87
|
1052
|
Chris@87
|
1053 """
|
Chris@87
|
1054 c = lagval(x, c)
|
Chris@87
|
1055 c = lagval(y, c)
|
Chris@87
|
1056 return c
|
Chris@87
|
1057
|
Chris@87
|
1058
|
Chris@87
|
1059 def lagval3d(x, y, z, c):
|
Chris@87
|
1060 """
|
Chris@87
|
1061 Evaluate a 3-D Laguerre series at points (x, y, z).
|
Chris@87
|
1062
|
Chris@87
|
1063 This function returns the values:
|
Chris@87
|
1064
|
Chris@87
|
1065 .. math:: p(x,y,z) = \\sum_{i,j,k} c_{i,j,k} * L_i(x) * L_j(y) * L_k(z)
|
Chris@87
|
1066
|
Chris@87
|
1067 The parameters `x`, `y`, and `z` are converted to arrays only if
|
Chris@87
|
1068 they are tuples or a lists, otherwise they are treated as a scalars and
|
Chris@87
|
1069 they must have the same shape after conversion. In either case, either
|
Chris@87
|
1070 `x`, `y`, and `z` or their elements must support multiplication and
|
Chris@87
|
1071 addition both with themselves and with the elements of `c`.
|
Chris@87
|
1072
|
Chris@87
|
1073 If `c` has fewer than 3 dimensions, ones are implicitly appended to its
|
Chris@87
|
1074 shape to make it 3-D. The shape of the result will be c.shape[3:] +
|
Chris@87
|
1075 x.shape.
|
Chris@87
|
1076
|
Chris@87
|
1077 Parameters
|
Chris@87
|
1078 ----------
|
Chris@87
|
1079 x, y, z : array_like, compatible object
|
Chris@87
|
1080 The three dimensional series is evaluated at the points
|
Chris@87
|
1081 `(x, y, z)`, where `x`, `y`, and `z` must have the same shape. If
|
Chris@87
|
1082 any of `x`, `y`, or `z` is a list or tuple, it is first converted
|
Chris@87
|
1083 to an ndarray, otherwise it is left unchanged and if it isn't an
|
Chris@87
|
1084 ndarray it is treated as a scalar.
|
Chris@87
|
1085 c : array_like
|
Chris@87
|
1086 Array of coefficients ordered so that the coefficient of the term of
|
Chris@87
|
1087 multi-degree i,j,k is contained in ``c[i,j,k]``. If `c` has dimension
|
Chris@87
|
1088 greater than 3 the remaining indices enumerate multiple sets of
|
Chris@87
|
1089 coefficients.
|
Chris@87
|
1090
|
Chris@87
|
1091 Returns
|
Chris@87
|
1092 -------
|
Chris@87
|
1093 values : ndarray, compatible object
|
Chris@87
|
1094 The values of the multidimension polynomial on points formed with
|
Chris@87
|
1095 triples of corresponding values from `x`, `y`, and `z`.
|
Chris@87
|
1096
|
Chris@87
|
1097 See Also
|
Chris@87
|
1098 --------
|
Chris@87
|
1099 lagval, lagval2d, laggrid2d, laggrid3d
|
Chris@87
|
1100
|
Chris@87
|
1101 Notes
|
Chris@87
|
1102 -----
|
Chris@87
|
1103
|
Chris@87
|
1104 .. versionadded::1.7.0
|
Chris@87
|
1105
|
Chris@87
|
1106 """
|
Chris@87
|
1107 try:
|
Chris@87
|
1108 x, y, z = np.array((x, y, z), copy=0)
|
Chris@87
|
1109 except:
|
Chris@87
|
1110 raise ValueError('x, y, z are incompatible')
|
Chris@87
|
1111
|
Chris@87
|
1112 c = lagval(x, c)
|
Chris@87
|
1113 c = lagval(y, c, tensor=False)
|
Chris@87
|
1114 c = lagval(z, c, tensor=False)
|
Chris@87
|
1115 return c
|
Chris@87
|
1116
|
Chris@87
|
1117
|
Chris@87
|
1118 def laggrid3d(x, y, z, c):
|
Chris@87
|
1119 """
|
Chris@87
|
1120 Evaluate a 3-D Laguerre series on the Cartesian product of x, y, and z.
|
Chris@87
|
1121
|
Chris@87
|
1122 This function returns the values:
|
Chris@87
|
1123
|
Chris@87
|
1124 .. math:: p(a,b,c) = \\sum_{i,j,k} c_{i,j,k} * L_i(a) * L_j(b) * L_k(c)
|
Chris@87
|
1125
|
Chris@87
|
1126 where the points `(a, b, c)` consist of all triples formed by taking
|
Chris@87
|
1127 `a` from `x`, `b` from `y`, and `c` from `z`. The resulting points form
|
Chris@87
|
1128 a grid with `x` in the first dimension, `y` in the second, and `z` in
|
Chris@87
|
1129 the third.
|
Chris@87
|
1130
|
Chris@87
|
1131 The parameters `x`, `y`, and `z` are converted to arrays only if they
|
Chris@87
|
1132 are tuples or a lists, otherwise they are treated as a scalars. In
|
Chris@87
|
1133 either case, either `x`, `y`, and `z` or their elements must support
|
Chris@87
|
1134 multiplication and addition both with themselves and with the elements
|
Chris@87
|
1135 of `c`.
|
Chris@87
|
1136
|
Chris@87
|
1137 If `c` has fewer than three dimensions, ones are implicitly appended to
|
Chris@87
|
1138 its shape to make it 3-D. The shape of the result will be c.shape[3:] +
|
Chris@87
|
1139 x.shape + y.shape + z.shape.
|
Chris@87
|
1140
|
Chris@87
|
1141 Parameters
|
Chris@87
|
1142 ----------
|
Chris@87
|
1143 x, y, z : array_like, compatible objects
|
Chris@87
|
1144 The three dimensional series is evaluated at the points in the
|
Chris@87
|
1145 Cartesian product of `x`, `y`, and `z`. If `x`,`y`, or `z` is a
|
Chris@87
|
1146 list or tuple, it is first converted to an ndarray, otherwise it is
|
Chris@87
|
1147 left unchanged and, if it isn't an ndarray, it is treated as a
|
Chris@87
|
1148 scalar.
|
Chris@87
|
1149 c : array_like
|
Chris@87
|
1150 Array of coefficients ordered so that the coefficients for terms of
|
Chris@87
|
1151 degree i,j are contained in ``c[i,j]``. If `c` has dimension
|
Chris@87
|
1152 greater than two the remaining indices enumerate multiple sets of
|
Chris@87
|
1153 coefficients.
|
Chris@87
|
1154
|
Chris@87
|
1155 Returns
|
Chris@87
|
1156 -------
|
Chris@87
|
1157 values : ndarray, compatible object
|
Chris@87
|
1158 The values of the two dimensional polynomial at points in the Cartesian
|
Chris@87
|
1159 product of `x` and `y`.
|
Chris@87
|
1160
|
Chris@87
|
1161 See Also
|
Chris@87
|
1162 --------
|
Chris@87
|
1163 lagval, lagval2d, laggrid2d, lagval3d
|
Chris@87
|
1164
|
Chris@87
|
1165 Notes
|
Chris@87
|
1166 -----
|
Chris@87
|
1167
|
Chris@87
|
1168 .. versionadded::1.7.0
|
Chris@87
|
1169
|
Chris@87
|
1170 """
|
Chris@87
|
1171 c = lagval(x, c)
|
Chris@87
|
1172 c = lagval(y, c)
|
Chris@87
|
1173 c = lagval(z, c)
|
Chris@87
|
1174 return c
|
Chris@87
|
1175
|
Chris@87
|
1176
|
Chris@87
|
1177 def lagvander(x, deg):
|
Chris@87
|
1178 """Pseudo-Vandermonde matrix of given degree.
|
Chris@87
|
1179
|
Chris@87
|
1180 Returns the pseudo-Vandermonde matrix of degree `deg` and sample points
|
Chris@87
|
1181 `x`. The pseudo-Vandermonde matrix is defined by
|
Chris@87
|
1182
|
Chris@87
|
1183 .. math:: V[..., i] = L_i(x)
|
Chris@87
|
1184
|
Chris@87
|
1185 where `0 <= i <= deg`. The leading indices of `V` index the elements of
|
Chris@87
|
1186 `x` and the last index is the degree of the Laguerre polynomial.
|
Chris@87
|
1187
|
Chris@87
|
1188 If `c` is a 1-D array of coefficients of length `n + 1` and `V` is the
|
Chris@87
|
1189 array ``V = lagvander(x, n)``, then ``np.dot(V, c)`` and
|
Chris@87
|
1190 ``lagval(x, c)`` are the same up to roundoff. This equivalence is
|
Chris@87
|
1191 useful both for least squares fitting and for the evaluation of a large
|
Chris@87
|
1192 number of Laguerre series of the same degree and sample points.
|
Chris@87
|
1193
|
Chris@87
|
1194 Parameters
|
Chris@87
|
1195 ----------
|
Chris@87
|
1196 x : array_like
|
Chris@87
|
1197 Array of points. The dtype is converted to float64 or complex128
|
Chris@87
|
1198 depending on whether any of the elements are complex. If `x` is
|
Chris@87
|
1199 scalar it is converted to a 1-D array.
|
Chris@87
|
1200 deg : int
|
Chris@87
|
1201 Degree of the resulting matrix.
|
Chris@87
|
1202
|
Chris@87
|
1203 Returns
|
Chris@87
|
1204 -------
|
Chris@87
|
1205 vander : ndarray
|
Chris@87
|
1206 The pseudo-Vandermonde matrix. The shape of the returned matrix is
|
Chris@87
|
1207 ``x.shape + (deg + 1,)``, where The last index is the degree of the
|
Chris@87
|
1208 corresponding Laguerre polynomial. The dtype will be the same as
|
Chris@87
|
1209 the converted `x`.
|
Chris@87
|
1210
|
Chris@87
|
1211 Examples
|
Chris@87
|
1212 --------
|
Chris@87
|
1213 >>> from numpy.polynomial.laguerre import lagvander
|
Chris@87
|
1214 >>> x = np.array([0, 1, 2])
|
Chris@87
|
1215 >>> lagvander(x, 3)
|
Chris@87
|
1216 array([[ 1. , 1. , 1. , 1. ],
|
Chris@87
|
1217 [ 1. , 0. , -0.5 , -0.66666667],
|
Chris@87
|
1218 [ 1. , -1. , -1. , -0.33333333]])
|
Chris@87
|
1219
|
Chris@87
|
1220 """
|
Chris@87
|
1221 ideg = int(deg)
|
Chris@87
|
1222 if ideg != deg:
|
Chris@87
|
1223 raise ValueError("deg must be integer")
|
Chris@87
|
1224 if ideg < 0:
|
Chris@87
|
1225 raise ValueError("deg must be non-negative")
|
Chris@87
|
1226
|
Chris@87
|
1227 x = np.array(x, copy=0, ndmin=1) + 0.0
|
Chris@87
|
1228 dims = (ideg + 1,) + x.shape
|
Chris@87
|
1229 dtyp = x.dtype
|
Chris@87
|
1230 v = np.empty(dims, dtype=dtyp)
|
Chris@87
|
1231 v[0] = x*0 + 1
|
Chris@87
|
1232 if ideg > 0:
|
Chris@87
|
1233 v[1] = 1 - x
|
Chris@87
|
1234 for i in range(2, ideg + 1):
|
Chris@87
|
1235 v[i] = (v[i-1]*(2*i - 1 - x) - v[i-2]*(i - 1))/i
|
Chris@87
|
1236 return np.rollaxis(v, 0, v.ndim)
|
Chris@87
|
1237
|
Chris@87
|
1238
|
Chris@87
|
1239 def lagvander2d(x, y, deg):
|
Chris@87
|
1240 """Pseudo-Vandermonde matrix of given degrees.
|
Chris@87
|
1241
|
Chris@87
|
1242 Returns the pseudo-Vandermonde matrix of degrees `deg` and sample
|
Chris@87
|
1243 points `(x, y)`. The pseudo-Vandermonde matrix is defined by
|
Chris@87
|
1244
|
Chris@87
|
1245 .. math:: V[..., deg[1]*i + j] = L_i(x) * L_j(y),
|
Chris@87
|
1246
|
Chris@87
|
1247 where `0 <= i <= deg[0]` and `0 <= j <= deg[1]`. The leading indices of
|
Chris@87
|
1248 `V` index the points `(x, y)` and the last index encodes the degrees of
|
Chris@87
|
1249 the Laguerre polynomials.
|
Chris@87
|
1250
|
Chris@87
|
1251 If ``V = lagvander2d(x, y, [xdeg, ydeg])``, then the columns of `V`
|
Chris@87
|
1252 correspond to the elements of a 2-D coefficient array `c` of shape
|
Chris@87
|
1253 (xdeg + 1, ydeg + 1) in the order
|
Chris@87
|
1254
|
Chris@87
|
1255 .. math:: c_{00}, c_{01}, c_{02} ... , c_{10}, c_{11}, c_{12} ...
|
Chris@87
|
1256
|
Chris@87
|
1257 and ``np.dot(V, c.flat)`` and ``lagval2d(x, y, c)`` will be the same
|
Chris@87
|
1258 up to roundoff. This equivalence is useful both for least squares
|
Chris@87
|
1259 fitting and for the evaluation of a large number of 2-D Laguerre
|
Chris@87
|
1260 series of the same degrees and sample points.
|
Chris@87
|
1261
|
Chris@87
|
1262 Parameters
|
Chris@87
|
1263 ----------
|
Chris@87
|
1264 x, y : array_like
|
Chris@87
|
1265 Arrays of point coordinates, all of the same shape. The dtypes
|
Chris@87
|
1266 will be converted to either float64 or complex128 depending on
|
Chris@87
|
1267 whether any of the elements are complex. Scalars are converted to
|
Chris@87
|
1268 1-D arrays.
|
Chris@87
|
1269 deg : list of ints
|
Chris@87
|
1270 List of maximum degrees of the form [x_deg, y_deg].
|
Chris@87
|
1271
|
Chris@87
|
1272 Returns
|
Chris@87
|
1273 -------
|
Chris@87
|
1274 vander2d : ndarray
|
Chris@87
|
1275 The shape of the returned matrix is ``x.shape + (order,)``, where
|
Chris@87
|
1276 :math:`order = (deg[0]+1)*(deg([1]+1)`. The dtype will be the same
|
Chris@87
|
1277 as the converted `x` and `y`.
|
Chris@87
|
1278
|
Chris@87
|
1279 See Also
|
Chris@87
|
1280 --------
|
Chris@87
|
1281 lagvander, lagvander3d. lagval2d, lagval3d
|
Chris@87
|
1282
|
Chris@87
|
1283 Notes
|
Chris@87
|
1284 -----
|
Chris@87
|
1285
|
Chris@87
|
1286 .. versionadded::1.7.0
|
Chris@87
|
1287
|
Chris@87
|
1288 """
|
Chris@87
|
1289 ideg = [int(d) for d in deg]
|
Chris@87
|
1290 is_valid = [id == d and id >= 0 for id, d in zip(ideg, deg)]
|
Chris@87
|
1291 if is_valid != [1, 1]:
|
Chris@87
|
1292 raise ValueError("degrees must be non-negative integers")
|
Chris@87
|
1293 degx, degy = ideg
|
Chris@87
|
1294 x, y = np.array((x, y), copy=0) + 0.0
|
Chris@87
|
1295
|
Chris@87
|
1296 vx = lagvander(x, degx)
|
Chris@87
|
1297 vy = lagvander(y, degy)
|
Chris@87
|
1298 v = vx[..., None]*vy[..., None,:]
|
Chris@87
|
1299 return v.reshape(v.shape[:-2] + (-1,))
|
Chris@87
|
1300
|
Chris@87
|
1301
|
Chris@87
|
1302 def lagvander3d(x, y, z, deg):
|
Chris@87
|
1303 """Pseudo-Vandermonde matrix of given degrees.
|
Chris@87
|
1304
|
Chris@87
|
1305 Returns the pseudo-Vandermonde matrix of degrees `deg` and sample
|
Chris@87
|
1306 points `(x, y, z)`. If `l, m, n` are the given degrees in `x, y, z`,
|
Chris@87
|
1307 then The pseudo-Vandermonde matrix is defined by
|
Chris@87
|
1308
|
Chris@87
|
1309 .. math:: V[..., (m+1)(n+1)i + (n+1)j + k] = L_i(x)*L_j(y)*L_k(z),
|
Chris@87
|
1310
|
Chris@87
|
1311 where `0 <= i <= l`, `0 <= j <= m`, and `0 <= j <= n`. The leading
|
Chris@87
|
1312 indices of `V` index the points `(x, y, z)` and the last index encodes
|
Chris@87
|
1313 the degrees of the Laguerre polynomials.
|
Chris@87
|
1314
|
Chris@87
|
1315 If ``V = lagvander3d(x, y, z, [xdeg, ydeg, zdeg])``, then the columns
|
Chris@87
|
1316 of `V` correspond to the elements of a 3-D coefficient array `c` of
|
Chris@87
|
1317 shape (xdeg + 1, ydeg + 1, zdeg + 1) in the order
|
Chris@87
|
1318
|
Chris@87
|
1319 .. math:: c_{000}, c_{001}, c_{002},... , c_{010}, c_{011}, c_{012},...
|
Chris@87
|
1320
|
Chris@87
|
1321 and ``np.dot(V, c.flat)`` and ``lagval3d(x, y, z, c)`` will be the
|
Chris@87
|
1322 same up to roundoff. This equivalence is useful both for least squares
|
Chris@87
|
1323 fitting and for the evaluation of a large number of 3-D Laguerre
|
Chris@87
|
1324 series of the same degrees and sample points.
|
Chris@87
|
1325
|
Chris@87
|
1326 Parameters
|
Chris@87
|
1327 ----------
|
Chris@87
|
1328 x, y, z : array_like
|
Chris@87
|
1329 Arrays of point coordinates, all of the same shape. The dtypes will
|
Chris@87
|
1330 be converted to either float64 or complex128 depending on whether
|
Chris@87
|
1331 any of the elements are complex. Scalars are converted to 1-D
|
Chris@87
|
1332 arrays.
|
Chris@87
|
1333 deg : list of ints
|
Chris@87
|
1334 List of maximum degrees of the form [x_deg, y_deg, z_deg].
|
Chris@87
|
1335
|
Chris@87
|
1336 Returns
|
Chris@87
|
1337 -------
|
Chris@87
|
1338 vander3d : ndarray
|
Chris@87
|
1339 The shape of the returned matrix is ``x.shape + (order,)``, where
|
Chris@87
|
1340 :math:`order = (deg[0]+1)*(deg([1]+1)*(deg[2]+1)`. The dtype will
|
Chris@87
|
1341 be the same as the converted `x`, `y`, and `z`.
|
Chris@87
|
1342
|
Chris@87
|
1343 See Also
|
Chris@87
|
1344 --------
|
Chris@87
|
1345 lagvander, lagvander3d. lagval2d, lagval3d
|
Chris@87
|
1346
|
Chris@87
|
1347 Notes
|
Chris@87
|
1348 -----
|
Chris@87
|
1349
|
Chris@87
|
1350 .. versionadded::1.7.0
|
Chris@87
|
1351
|
Chris@87
|
1352 """
|
Chris@87
|
1353 ideg = [int(d) for d in deg]
|
Chris@87
|
1354 is_valid = [id == d and id >= 0 for id, d in zip(ideg, deg)]
|
Chris@87
|
1355 if is_valid != [1, 1, 1]:
|
Chris@87
|
1356 raise ValueError("degrees must be non-negative integers")
|
Chris@87
|
1357 degx, degy, degz = ideg
|
Chris@87
|
1358 x, y, z = np.array((x, y, z), copy=0) + 0.0
|
Chris@87
|
1359
|
Chris@87
|
1360 vx = lagvander(x, degx)
|
Chris@87
|
1361 vy = lagvander(y, degy)
|
Chris@87
|
1362 vz = lagvander(z, degz)
|
Chris@87
|
1363 v = vx[..., None, None]*vy[..., None,:, None]*vz[..., None, None,:]
|
Chris@87
|
1364 return v.reshape(v.shape[:-3] + (-1,))
|
Chris@87
|
1365
|
Chris@87
|
1366
|
Chris@87
|
1367 def lagfit(x, y, deg, rcond=None, full=False, w=None):
|
Chris@87
|
1368 """
|
Chris@87
|
1369 Least squares fit of Laguerre series to data.
|
Chris@87
|
1370
|
Chris@87
|
1371 Return the coefficients of a Laguerre series of degree `deg` that is the
|
Chris@87
|
1372 least squares fit to the data values `y` given at points `x`. If `y` is
|
Chris@87
|
1373 1-D the returned coefficients will also be 1-D. If `y` is 2-D multiple
|
Chris@87
|
1374 fits are done, one for each column of `y`, and the resulting
|
Chris@87
|
1375 coefficients are stored in the corresponding columns of a 2-D return.
|
Chris@87
|
1376 The fitted polynomial(s) are in the form
|
Chris@87
|
1377
|
Chris@87
|
1378 .. math:: p(x) = c_0 + c_1 * L_1(x) + ... + c_n * L_n(x),
|
Chris@87
|
1379
|
Chris@87
|
1380 where `n` is `deg`.
|
Chris@87
|
1381
|
Chris@87
|
1382 Parameters
|
Chris@87
|
1383 ----------
|
Chris@87
|
1384 x : array_like, shape (M,)
|
Chris@87
|
1385 x-coordinates of the M sample points ``(x[i], y[i])``.
|
Chris@87
|
1386 y : array_like, shape (M,) or (M, K)
|
Chris@87
|
1387 y-coordinates of the sample points. Several data sets of sample
|
Chris@87
|
1388 points sharing the same x-coordinates can be fitted at once by
|
Chris@87
|
1389 passing in a 2D-array that contains one dataset per column.
|
Chris@87
|
1390 deg : int
|
Chris@87
|
1391 Degree of the fitting polynomial
|
Chris@87
|
1392 rcond : float, optional
|
Chris@87
|
1393 Relative condition number of the fit. Singular values smaller than
|
Chris@87
|
1394 this relative to the largest singular value will be ignored. The
|
Chris@87
|
1395 default value is len(x)*eps, where eps is the relative precision of
|
Chris@87
|
1396 the float type, about 2e-16 in most cases.
|
Chris@87
|
1397 full : bool, optional
|
Chris@87
|
1398 Switch determining nature of return value. When it is False (the
|
Chris@87
|
1399 default) just the coefficients are returned, when True diagnostic
|
Chris@87
|
1400 information from the singular value decomposition is also returned.
|
Chris@87
|
1401 w : array_like, shape (`M`,), optional
|
Chris@87
|
1402 Weights. If not None, the contribution of each point
|
Chris@87
|
1403 ``(x[i],y[i])`` to the fit is weighted by `w[i]`. Ideally the
|
Chris@87
|
1404 weights are chosen so that the errors of the products ``w[i]*y[i]``
|
Chris@87
|
1405 all have the same variance. The default value is None.
|
Chris@87
|
1406
|
Chris@87
|
1407 Returns
|
Chris@87
|
1408 -------
|
Chris@87
|
1409 coef : ndarray, shape (M,) or (M, K)
|
Chris@87
|
1410 Laguerre coefficients ordered from low to high. If `y` was 2-D,
|
Chris@87
|
1411 the coefficients for the data in column k of `y` are in column
|
Chris@87
|
1412 `k`.
|
Chris@87
|
1413
|
Chris@87
|
1414 [residuals, rank, singular_values, rcond] : list
|
Chris@87
|
1415 These values are only returned if `full` = True
|
Chris@87
|
1416
|
Chris@87
|
1417 resid -- sum of squared residuals of the least squares fit
|
Chris@87
|
1418 rank -- the numerical rank of the scaled Vandermonde matrix
|
Chris@87
|
1419 sv -- singular values of the scaled Vandermonde matrix
|
Chris@87
|
1420 rcond -- value of `rcond`.
|
Chris@87
|
1421
|
Chris@87
|
1422 For more details, see `linalg.lstsq`.
|
Chris@87
|
1423
|
Chris@87
|
1424 Warns
|
Chris@87
|
1425 -----
|
Chris@87
|
1426 RankWarning
|
Chris@87
|
1427 The rank of the coefficient matrix in the least-squares fit is
|
Chris@87
|
1428 deficient. The warning is only raised if `full` = False. The
|
Chris@87
|
1429 warnings can be turned off by
|
Chris@87
|
1430
|
Chris@87
|
1431 >>> import warnings
|
Chris@87
|
1432 >>> warnings.simplefilter('ignore', RankWarning)
|
Chris@87
|
1433
|
Chris@87
|
1434 See Also
|
Chris@87
|
1435 --------
|
Chris@87
|
1436 chebfit, legfit, polyfit, hermfit, hermefit
|
Chris@87
|
1437 lagval : Evaluates a Laguerre series.
|
Chris@87
|
1438 lagvander : pseudo Vandermonde matrix of Laguerre series.
|
Chris@87
|
1439 lagweight : Laguerre weight function.
|
Chris@87
|
1440 linalg.lstsq : Computes a least-squares fit from the matrix.
|
Chris@87
|
1441 scipy.interpolate.UnivariateSpline : Computes spline fits.
|
Chris@87
|
1442
|
Chris@87
|
1443 Notes
|
Chris@87
|
1444 -----
|
Chris@87
|
1445 The solution is the coefficients of the Laguerre series `p` that
|
Chris@87
|
1446 minimizes the sum of the weighted squared errors
|
Chris@87
|
1447
|
Chris@87
|
1448 .. math:: E = \\sum_j w_j^2 * |y_j - p(x_j)|^2,
|
Chris@87
|
1449
|
Chris@87
|
1450 where the :math:`w_j` are the weights. This problem is solved by
|
Chris@87
|
1451 setting up as the (typically) overdetermined matrix equation
|
Chris@87
|
1452
|
Chris@87
|
1453 .. math:: V(x) * c = w * y,
|
Chris@87
|
1454
|
Chris@87
|
1455 where `V` is the weighted pseudo Vandermonde matrix of `x`, `c` are the
|
Chris@87
|
1456 coefficients to be solved for, `w` are the weights, and `y` are the
|
Chris@87
|
1457 observed values. This equation is then solved using the singular value
|
Chris@87
|
1458 decomposition of `V`.
|
Chris@87
|
1459
|
Chris@87
|
1460 If some of the singular values of `V` are so small that they are
|
Chris@87
|
1461 neglected, then a `RankWarning` will be issued. This means that the
|
Chris@87
|
1462 coefficient values may be poorly determined. Using a lower order fit
|
Chris@87
|
1463 will usually get rid of the warning. The `rcond` parameter can also be
|
Chris@87
|
1464 set to a value smaller than its default, but the resulting fit may be
|
Chris@87
|
1465 spurious and have large contributions from roundoff error.
|
Chris@87
|
1466
|
Chris@87
|
1467 Fits using Laguerre series are probably most useful when the data can
|
Chris@87
|
1468 be approximated by ``sqrt(w(x)) * p(x)``, where `w(x)` is the Laguerre
|
Chris@87
|
1469 weight. In that case the weight ``sqrt(w(x[i])`` should be used
|
Chris@87
|
1470 together with data values ``y[i]/sqrt(w(x[i])``. The weight function is
|
Chris@87
|
1471 available as `lagweight`.
|
Chris@87
|
1472
|
Chris@87
|
1473 References
|
Chris@87
|
1474 ----------
|
Chris@87
|
1475 .. [1] Wikipedia, "Curve fitting",
|
Chris@87
|
1476 http://en.wikipedia.org/wiki/Curve_fitting
|
Chris@87
|
1477
|
Chris@87
|
1478 Examples
|
Chris@87
|
1479 --------
|
Chris@87
|
1480 >>> from numpy.polynomial.laguerre import lagfit, lagval
|
Chris@87
|
1481 >>> x = np.linspace(0, 10)
|
Chris@87
|
1482 >>> err = np.random.randn(len(x))/10
|
Chris@87
|
1483 >>> y = lagval(x, [1, 2, 3]) + err
|
Chris@87
|
1484 >>> lagfit(x, y, 2)
|
Chris@87
|
1485 array([ 0.96971004, 2.00193749, 3.00288744])
|
Chris@87
|
1486
|
Chris@87
|
1487 """
|
Chris@87
|
1488 order = int(deg) + 1
|
Chris@87
|
1489 x = np.asarray(x) + 0.0
|
Chris@87
|
1490 y = np.asarray(y) + 0.0
|
Chris@87
|
1491
|
Chris@87
|
1492 # check arguments.
|
Chris@87
|
1493 if deg < 0:
|
Chris@87
|
1494 raise ValueError("expected deg >= 0")
|
Chris@87
|
1495 if x.ndim != 1:
|
Chris@87
|
1496 raise TypeError("expected 1D vector for x")
|
Chris@87
|
1497 if x.size == 0:
|
Chris@87
|
1498 raise TypeError("expected non-empty vector for x")
|
Chris@87
|
1499 if y.ndim < 1 or y.ndim > 2:
|
Chris@87
|
1500 raise TypeError("expected 1D or 2D array for y")
|
Chris@87
|
1501 if len(x) != len(y):
|
Chris@87
|
1502 raise TypeError("expected x and y to have same length")
|
Chris@87
|
1503
|
Chris@87
|
1504 # set up the least squares matrices in transposed form
|
Chris@87
|
1505 lhs = lagvander(x, deg).T
|
Chris@87
|
1506 rhs = y.T
|
Chris@87
|
1507 if w is not None:
|
Chris@87
|
1508 w = np.asarray(w) + 0.0
|
Chris@87
|
1509 if w.ndim != 1:
|
Chris@87
|
1510 raise TypeError("expected 1D vector for w")
|
Chris@87
|
1511 if len(x) != len(w):
|
Chris@87
|
1512 raise TypeError("expected x and w to have same length")
|
Chris@87
|
1513 # apply weights. Don't use inplace operations as they
|
Chris@87
|
1514 # can cause problems with NA.
|
Chris@87
|
1515 lhs = lhs * w
|
Chris@87
|
1516 rhs = rhs * w
|
Chris@87
|
1517
|
Chris@87
|
1518 # set rcond
|
Chris@87
|
1519 if rcond is None:
|
Chris@87
|
1520 rcond = len(x)*np.finfo(x.dtype).eps
|
Chris@87
|
1521
|
Chris@87
|
1522 # Determine the norms of the design matrix columns.
|
Chris@87
|
1523 if issubclass(lhs.dtype.type, np.complexfloating):
|
Chris@87
|
1524 scl = np.sqrt((np.square(lhs.real) + np.square(lhs.imag)).sum(1))
|
Chris@87
|
1525 else:
|
Chris@87
|
1526 scl = np.sqrt(np.square(lhs).sum(1))
|
Chris@87
|
1527 scl[scl == 0] = 1
|
Chris@87
|
1528
|
Chris@87
|
1529 # Solve the least squares problem.
|
Chris@87
|
1530 c, resids, rank, s = la.lstsq(lhs.T/scl, rhs.T, rcond)
|
Chris@87
|
1531 c = (c.T/scl).T
|
Chris@87
|
1532
|
Chris@87
|
1533 # warn on rank reduction
|
Chris@87
|
1534 if rank != order and not full:
|
Chris@87
|
1535 msg = "The fit may be poorly conditioned"
|
Chris@87
|
1536 warnings.warn(msg, pu.RankWarning)
|
Chris@87
|
1537
|
Chris@87
|
1538 if full:
|
Chris@87
|
1539 return c, [resids, rank, s, rcond]
|
Chris@87
|
1540 else:
|
Chris@87
|
1541 return c
|
Chris@87
|
1542
|
Chris@87
|
1543
|
Chris@87
|
1544 def lagcompanion(c):
|
Chris@87
|
1545 """
|
Chris@87
|
1546 Return the companion matrix of c.
|
Chris@87
|
1547
|
Chris@87
|
1548 The usual companion matrix of the Laguerre polynomials is already
|
Chris@87
|
1549 symmetric when `c` is a basis Laguerre polynomial, so no scaling is
|
Chris@87
|
1550 applied.
|
Chris@87
|
1551
|
Chris@87
|
1552 Parameters
|
Chris@87
|
1553 ----------
|
Chris@87
|
1554 c : array_like
|
Chris@87
|
1555 1-D array of Laguerre series coefficients ordered from low to high
|
Chris@87
|
1556 degree.
|
Chris@87
|
1557
|
Chris@87
|
1558 Returns
|
Chris@87
|
1559 -------
|
Chris@87
|
1560 mat : ndarray
|
Chris@87
|
1561 Companion matrix of dimensions (deg, deg).
|
Chris@87
|
1562
|
Chris@87
|
1563 Notes
|
Chris@87
|
1564 -----
|
Chris@87
|
1565
|
Chris@87
|
1566 .. versionadded::1.7.0
|
Chris@87
|
1567
|
Chris@87
|
1568 """
|
Chris@87
|
1569 # c is a trimmed copy
|
Chris@87
|
1570 [c] = pu.as_series([c])
|
Chris@87
|
1571 if len(c) < 2:
|
Chris@87
|
1572 raise ValueError('Series must have maximum degree of at least 1.')
|
Chris@87
|
1573 if len(c) == 2:
|
Chris@87
|
1574 return np.array([[1 + c[0]/c[1]]])
|
Chris@87
|
1575
|
Chris@87
|
1576 n = len(c) - 1
|
Chris@87
|
1577 mat = np.zeros((n, n), dtype=c.dtype)
|
Chris@87
|
1578 top = mat.reshape(-1)[1::n+1]
|
Chris@87
|
1579 mid = mat.reshape(-1)[0::n+1]
|
Chris@87
|
1580 bot = mat.reshape(-1)[n::n+1]
|
Chris@87
|
1581 top[...] = -np.arange(1, n)
|
Chris@87
|
1582 mid[...] = 2.*np.arange(n) + 1.
|
Chris@87
|
1583 bot[...] = top
|
Chris@87
|
1584 mat[:, -1] += (c[:-1]/c[-1])*n
|
Chris@87
|
1585 return mat
|
Chris@87
|
1586
|
Chris@87
|
1587
|
Chris@87
|
1588 def lagroots(c):
|
Chris@87
|
1589 """
|
Chris@87
|
1590 Compute the roots of a Laguerre series.
|
Chris@87
|
1591
|
Chris@87
|
1592 Return the roots (a.k.a. "zeros") of the polynomial
|
Chris@87
|
1593
|
Chris@87
|
1594 .. math:: p(x) = \\sum_i c[i] * L_i(x).
|
Chris@87
|
1595
|
Chris@87
|
1596 Parameters
|
Chris@87
|
1597 ----------
|
Chris@87
|
1598 c : 1-D array_like
|
Chris@87
|
1599 1-D array of coefficients.
|
Chris@87
|
1600
|
Chris@87
|
1601 Returns
|
Chris@87
|
1602 -------
|
Chris@87
|
1603 out : ndarray
|
Chris@87
|
1604 Array of the roots of the series. If all the roots are real,
|
Chris@87
|
1605 then `out` is also real, otherwise it is complex.
|
Chris@87
|
1606
|
Chris@87
|
1607 See Also
|
Chris@87
|
1608 --------
|
Chris@87
|
1609 polyroots, legroots, chebroots, hermroots, hermeroots
|
Chris@87
|
1610
|
Chris@87
|
1611 Notes
|
Chris@87
|
1612 -----
|
Chris@87
|
1613 The root estimates are obtained as the eigenvalues of the companion
|
Chris@87
|
1614 matrix, Roots far from the origin of the complex plane may have large
|
Chris@87
|
1615 errors due to the numerical instability of the series for such
|
Chris@87
|
1616 values. Roots with multiplicity greater than 1 will also show larger
|
Chris@87
|
1617 errors as the value of the series near such points is relatively
|
Chris@87
|
1618 insensitive to errors in the roots. Isolated roots near the origin can
|
Chris@87
|
1619 be improved by a few iterations of Newton's method.
|
Chris@87
|
1620
|
Chris@87
|
1621 The Laguerre series basis polynomials aren't powers of `x` so the
|
Chris@87
|
1622 results of this function may seem unintuitive.
|
Chris@87
|
1623
|
Chris@87
|
1624 Examples
|
Chris@87
|
1625 --------
|
Chris@87
|
1626 >>> from numpy.polynomial.laguerre import lagroots, lagfromroots
|
Chris@87
|
1627 >>> coef = lagfromroots([0, 1, 2])
|
Chris@87
|
1628 >>> coef
|
Chris@87
|
1629 array([ 2., -8., 12., -6.])
|
Chris@87
|
1630 >>> lagroots(coef)
|
Chris@87
|
1631 array([ -4.44089210e-16, 1.00000000e+00, 2.00000000e+00])
|
Chris@87
|
1632
|
Chris@87
|
1633 """
|
Chris@87
|
1634 # c is a trimmed copy
|
Chris@87
|
1635 [c] = pu.as_series([c])
|
Chris@87
|
1636 if len(c) <= 1:
|
Chris@87
|
1637 return np.array([], dtype=c.dtype)
|
Chris@87
|
1638 if len(c) == 2:
|
Chris@87
|
1639 return np.array([1 + c[0]/c[1]])
|
Chris@87
|
1640
|
Chris@87
|
1641 m = lagcompanion(c)
|
Chris@87
|
1642 r = la.eigvals(m)
|
Chris@87
|
1643 r.sort()
|
Chris@87
|
1644 return r
|
Chris@87
|
1645
|
Chris@87
|
1646
|
Chris@87
|
1647 def laggauss(deg):
|
Chris@87
|
1648 """
|
Chris@87
|
1649 Gauss-Laguerre quadrature.
|
Chris@87
|
1650
|
Chris@87
|
1651 Computes the sample points and weights for Gauss-Laguerre quadrature.
|
Chris@87
|
1652 These sample points and weights will correctly integrate polynomials of
|
Chris@87
|
1653 degree :math:`2*deg - 1` or less over the interval :math:`[0, \inf]`
|
Chris@87
|
1654 with the weight function :math:`f(x) = \exp(-x)`.
|
Chris@87
|
1655
|
Chris@87
|
1656 Parameters
|
Chris@87
|
1657 ----------
|
Chris@87
|
1658 deg : int
|
Chris@87
|
1659 Number of sample points and weights. It must be >= 1.
|
Chris@87
|
1660
|
Chris@87
|
1661 Returns
|
Chris@87
|
1662 -------
|
Chris@87
|
1663 x : ndarray
|
Chris@87
|
1664 1-D ndarray containing the sample points.
|
Chris@87
|
1665 y : ndarray
|
Chris@87
|
1666 1-D ndarray containing the weights.
|
Chris@87
|
1667
|
Chris@87
|
1668 Notes
|
Chris@87
|
1669 -----
|
Chris@87
|
1670
|
Chris@87
|
1671 .. versionadded::1.7.0
|
Chris@87
|
1672
|
Chris@87
|
1673 The results have only been tested up to degree 100 higher degrees may
|
Chris@87
|
1674 be problematic. The weights are determined by using the fact that
|
Chris@87
|
1675
|
Chris@87
|
1676 .. math:: w_k = c / (L'_n(x_k) * L_{n-1}(x_k))
|
Chris@87
|
1677
|
Chris@87
|
1678 where :math:`c` is a constant independent of :math:`k` and :math:`x_k`
|
Chris@87
|
1679 is the k'th root of :math:`L_n`, and then scaling the results to get
|
Chris@87
|
1680 the right value when integrating 1.
|
Chris@87
|
1681
|
Chris@87
|
1682 """
|
Chris@87
|
1683 ideg = int(deg)
|
Chris@87
|
1684 if ideg != deg or ideg < 1:
|
Chris@87
|
1685 raise ValueError("deg must be a non-negative integer")
|
Chris@87
|
1686
|
Chris@87
|
1687 # first approximation of roots. We use the fact that the companion
|
Chris@87
|
1688 # matrix is symmetric in this case in order to obtain better zeros.
|
Chris@87
|
1689 c = np.array([0]*deg + [1])
|
Chris@87
|
1690 m = lagcompanion(c)
|
Chris@87
|
1691 x = la.eigvals(m)
|
Chris@87
|
1692 x.sort()
|
Chris@87
|
1693
|
Chris@87
|
1694 # improve roots by one application of Newton
|
Chris@87
|
1695 dy = lagval(x, c)
|
Chris@87
|
1696 df = lagval(x, lagder(c))
|
Chris@87
|
1697 x -= dy/df
|
Chris@87
|
1698
|
Chris@87
|
1699 # compute the weights. We scale the factor to avoid possible numerical
|
Chris@87
|
1700 # overflow.
|
Chris@87
|
1701 fm = lagval(x, c[1:])
|
Chris@87
|
1702 fm /= np.abs(fm).max()
|
Chris@87
|
1703 df /= np.abs(df).max()
|
Chris@87
|
1704 w = 1/(fm * df)
|
Chris@87
|
1705
|
Chris@87
|
1706 # scale w to get the right value, 1 in this case
|
Chris@87
|
1707 w /= w.sum()
|
Chris@87
|
1708
|
Chris@87
|
1709 return x, w
|
Chris@87
|
1710
|
Chris@87
|
1711
|
Chris@87
|
1712 def lagweight(x):
|
Chris@87
|
1713 """Weight function of the Laguerre polynomials.
|
Chris@87
|
1714
|
Chris@87
|
1715 The weight function is :math:`exp(-x)` and the interval of integration
|
Chris@87
|
1716 is :math:`[0, \inf]`. The Laguerre polynomials are orthogonal, but not
|
Chris@87
|
1717 normalized, with respect to this weight function.
|
Chris@87
|
1718
|
Chris@87
|
1719 Parameters
|
Chris@87
|
1720 ----------
|
Chris@87
|
1721 x : array_like
|
Chris@87
|
1722 Values at which the weight function will be computed.
|
Chris@87
|
1723
|
Chris@87
|
1724 Returns
|
Chris@87
|
1725 -------
|
Chris@87
|
1726 w : ndarray
|
Chris@87
|
1727 The weight function at `x`.
|
Chris@87
|
1728
|
Chris@87
|
1729 Notes
|
Chris@87
|
1730 -----
|
Chris@87
|
1731
|
Chris@87
|
1732 .. versionadded::1.7.0
|
Chris@87
|
1733
|
Chris@87
|
1734 """
|
Chris@87
|
1735 w = np.exp(-x)
|
Chris@87
|
1736 return w
|
Chris@87
|
1737
|
Chris@87
|
1738 #
|
Chris@87
|
1739 # Laguerre series class
|
Chris@87
|
1740 #
|
Chris@87
|
1741
|
Chris@87
|
1742 class Laguerre(ABCPolyBase):
|
Chris@87
|
1743 """A Laguerre series class.
|
Chris@87
|
1744
|
Chris@87
|
1745 The Laguerre class provides the standard Python numerical methods
|
Chris@87
|
1746 '+', '-', '*', '//', '%', 'divmod', '**', and '()' as well as the
|
Chris@87
|
1747 attributes and methods listed in the `ABCPolyBase` documentation.
|
Chris@87
|
1748
|
Chris@87
|
1749 Parameters
|
Chris@87
|
1750 ----------
|
Chris@87
|
1751 coef : array_like
|
Chris@87
|
1752 Laguerre coefficients in order of increasing degree, i.e,
|
Chris@87
|
1753 ``(1, 2, 3)`` gives ``1*L_0(x) + 2*L_1(X) + 3*L_2(x)``.
|
Chris@87
|
1754 domain : (2,) array_like, optional
|
Chris@87
|
1755 Domain to use. The interval ``[domain[0], domain[1]]`` is mapped
|
Chris@87
|
1756 to the interval ``[window[0], window[1]]`` by shifting and scaling.
|
Chris@87
|
1757 The default value is [0, 1].
|
Chris@87
|
1758 window : (2,) array_like, optional
|
Chris@87
|
1759 Window, see `domain` for its use. The default value is [0, 1].
|
Chris@87
|
1760
|
Chris@87
|
1761 .. versionadded:: 1.6.0
|
Chris@87
|
1762
|
Chris@87
|
1763 """
|
Chris@87
|
1764 # Virtual Functions
|
Chris@87
|
1765 _add = staticmethod(lagadd)
|
Chris@87
|
1766 _sub = staticmethod(lagsub)
|
Chris@87
|
1767 _mul = staticmethod(lagmul)
|
Chris@87
|
1768 _div = staticmethod(lagdiv)
|
Chris@87
|
1769 _pow = staticmethod(lagpow)
|
Chris@87
|
1770 _val = staticmethod(lagval)
|
Chris@87
|
1771 _int = staticmethod(lagint)
|
Chris@87
|
1772 _der = staticmethod(lagder)
|
Chris@87
|
1773 _fit = staticmethod(lagfit)
|
Chris@87
|
1774 _line = staticmethod(lagline)
|
Chris@87
|
1775 _roots = staticmethod(lagroots)
|
Chris@87
|
1776 _fromroots = staticmethod(lagfromroots)
|
Chris@87
|
1777
|
Chris@87
|
1778 # Virtual properties
|
Chris@87
|
1779 nickname = 'lag'
|
Chris@87
|
1780 domain = np.array(lagdomain)
|
Chris@87
|
1781 window = np.array(lagdomain)
|