Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/functions/piecewise.py
4057 views
1
r"""
2
Piecewise-defined Functions
3
4
Sage implements a very simple class of piecewise-defined functions.
5
Functions may be any type of symbolic expression. Infinite
6
intervals are not supported. The endpoints of each interval must
7
line up.
8
9
TODO:
10
11
- Implement max/min location and values,
12
13
- Need: parent object - ring of piecewise functions
14
15
- This class should derive from an element-type class, and should
16
define ``_add_``, ``_mul_``, etc. That will automatically take care
17
of left multiplication and proper coercion. The coercion mentioned
18
below for scalar mult on right is bad, since it only allows ints and
19
rationals. The right way is to use an element class and only define
20
``_mul_``, and have a parent, so anything gets coerced properly.
21
22
AUTHORS:
23
24
- David Joyner (2006-04): initial version
25
26
- David Joyner (2006-09): added __eq__, extend_by_zero_to, unextend,
27
convolution, trapezoid, trapezoid_integral_approximation,
28
riemann_sum, riemann_sum_integral_approximation, tangent_line fixed
29
bugs in __mul__, __add__
30
31
- David Joyner (2007-03): adding Hann filter for FS, added general FS
32
filter methods for computing and plotting, added options to plotting
33
of FS (eg, specifying rgb values are now allowed). Fixed bug in
34
documentation reported by Pablo De Napoli.
35
36
- David Joyner (2007-09): bug fixes due to behaviour of
37
SymbolicArithmetic
38
39
- David Joyner (2008-04): fixed docstring bugs reported by J Morrow; added
40
support for Laplace transform of functions with infinite support.
41
42
- David Joyner (2008-07): fixed a left multiplication bug reported by
43
C. Boncelet (by defining __rmul__ = __mul__).
44
45
- Paul Butler (2009-01): added indefinite integration and default_variable
46
47
TESTS::
48
49
sage: R.<x> = QQ[]
50
sage: f = Piecewise([[(0,1),1*x^0]])
51
sage: 2*f
52
Piecewise defined function with 1 parts, [[(0, 1), 2]]
53
"""
54
55
#*****************************************************************************
56
# Copyright (C) 2006 William Stein <[email protected]>
57
# 2006 David Joyner <[email protected]>
58
#
59
# Distributed under the terms of the GNU General Public License (GPL)
60
#
61
# This code is distributed in the hope that it will be useful,
62
# but WITHOUT ANY WARRANTY; without even the implied warranty of
63
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
64
# General Public License for more details.
65
#
66
# The full text of the GPL is available at:
67
#
68
# http://www.gnu.org/licenses/
69
#*****************************************************************************
70
71
from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing
72
from sage.misc.sage_eval import sage_eval
73
from sage.rings.all import QQ, RR, Integer, Rational, infinity
74
from sage.calculus.functional import derivative
75
from sage.symbolic.expression import is_Expression
76
from sage.symbolic.assumptions import assume, forget
77
78
from sage.calculus.calculus import SR, maxima
79
from sage.calculus.all import var
80
81
def piecewise(list_of_pairs, var=None):
82
"""
83
Returns a piecewise function from a list of (interval, function)
84
pairs.
85
86
``list_of_pairs`` is a list of pairs (I, fcn), where
87
fcn is a Sage function (such as a polynomial over RR, or functions
88
using the lambda notation), and I is an interval such as I = (1,3).
89
Two consecutive intervals must share a common endpoint.
90
91
If the optional ``var`` is specified, then any symbolic expressions
92
in the list will be converted to symbolic functions using
93
``fcn.function(var)``. (This says which variable is considered to
94
be "piecewise".)
95
96
We assume that these definitions are consistent (ie, no checking is
97
done).
98
99
EXAMPLES::
100
101
sage: f1(x) = -1
102
sage: f2(x) = 2
103
sage: f = Piecewise([[(0,pi/2),f1],[(pi/2,pi),f2]])
104
sage: f(1)
105
-1
106
sage: f(3)
107
2
108
sage: f = Piecewise([[(0,1),x], [(1,2),x^2]], x); f
109
Piecewise defined function with 2 parts, [[(0, 1), x |--> x], [(1, 2), x |--> x^2]]
110
sage: f(0.9)
111
0.900000000000000
112
sage: f(1.1)
113
1.21000000000000
114
"""
115
return PiecewisePolynomial(list_of_pairs, var=var)
116
117
Piecewise = piecewise
118
119
class PiecewisePolynomial:
120
"""
121
Returns a piecewise function from a list of (interval, function)
122
pairs.
123
124
EXAMPLES::
125
126
sage: f1(x) = -1
127
sage: f2(x) = 2
128
sage: f = Piecewise([[(0,pi/2),f1],[(pi/2,pi),f2]])
129
sage: f(1)
130
-1
131
sage: f(3)
132
2
133
"""
134
def __init__(self, list_of_pairs, var=None):
135
r"""
136
``list_of_pairs`` is a list of pairs (I, fcn), where
137
fcn is a Sage function (such as a polynomial over RR, or functions
138
using the lambda notation), and I is an interval such as I = (1,3).
139
Two consecutive intervals must share a common endpoint.
140
141
If the optional ``var`` is specified, then any symbolic
142
expressions in the list will be converted to symbolic
143
functions using ``fcn.function(var)``. (This says which
144
variable is considered to be "piecewise".)
145
146
We assume that these definitions are consistent (ie, no checking is
147
done).
148
149
EXAMPLES::
150
151
sage: f1(x) = 1
152
sage: f2(x) = 1 - x
153
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
154
sage: f.list()
155
[[(0, 1), x |--> 1], [(1, 2), x |--> -x + 1]]
156
sage: f.length()
157
2
158
"""
159
self._length = len(list_of_pairs)
160
self._intervals = [x[0] for x in list_of_pairs]
161
functions = [x[1] for x in list_of_pairs]
162
if var is not None:
163
for i in range(len(functions)):
164
if is_Expression(functions[i]):
165
functions[i] = functions[i].function(var)
166
self._functions = functions
167
# We regenerate self._list in case self._functions was modified
168
# above. This also protects us in case somebody mutates a list
169
# after they use it as an argument to piecewise().
170
self._list = [[self._intervals[i], self._functions[i]] for i in range(self._length)]
171
172
def list(self):
173
"""
174
Returns the pieces of this function as a list of functions.
175
176
EXAMPLES::
177
178
sage: f1(x) = 1
179
sage: f2(x) = 1 - x
180
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
181
sage: f.list()
182
[[(0, 1), x |--> 1], [(1, 2), x |--> -x + 1]]
183
"""
184
return self._list
185
186
def length(self):
187
"""
188
Returns the number of pieces of this function.
189
190
EXAMPLES::
191
192
sage: f1(x) = 1
193
sage: f2(x) = 1 - x
194
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
195
sage: f.length()
196
2
197
"""
198
return self._length
199
200
def __repr__(self):
201
"""
202
EXAMPLES::
203
204
sage: f1(x) = 1
205
sage: f2(x) = 1 - x
206
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]]); f
207
Piecewise defined function with 2 parts, [[(0, 1), x |--> 1], [(1, 2), x |--> -x + 1]]
208
"""
209
return 'Piecewise defined function with %s parts, %s'%(
210
self.length(),self.list())
211
212
def _latex_(self):
213
r"""
214
EXAMPLES::
215
216
sage: f1(x) = 1
217
sage: f2(x) = 1 - x
218
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
219
sage: latex(f)
220
\begin{cases}
221
x \ {\mapsto}\ 1 &\text{on $(0, 1)$}\cr
222
x \ {\mapsto}\ -x + 1 &\text{on $(1, 2)$}\cr
223
\end{cases}
224
225
::
226
227
sage: f(x) = sin(x*pi/2)
228
sage: g(x) = 1-(x-1)^2
229
sage: h(x) = -x
230
sage: P = Piecewise([[(0,1), f], [(1,3),g], [(3,5), h]])
231
sage: latex(P)
232
\begin{cases}
233
x \ {\mapsto}\ \sin\left(\frac{1}{2} \, \pi x\right) &\text{on $(0, 1)$}\cr
234
x \ {\mapsto}\ -{\left(x - 1\right)}^{2} + 1 &\text{on $(1, 3)$}\cr
235
x \ {\mapsto}\ -x &\text{on $(3, 5)$}\cr
236
\end{cases}
237
"""
238
from sage.misc.latex import latex
239
tex = ['\\begin{cases}\n']
240
for (left, right), f in self.list():
241
tex.append('%s &\\text{on $(%s, %s)$}\\cr\n' % (latex(f), left, right))
242
tex.append(r'\end{cases}')
243
return ''.join(tex)
244
245
def intervals(self):
246
"""
247
A piecewise non-polynomial example.
248
249
EXAMPLES::
250
251
sage: f1(x) = 1
252
sage: f2(x) = 1-x
253
sage: f3(x) = exp(x)
254
sage: f4(x) = sin(2*x)
255
sage: f = Piecewise([[(0,1),f1],[(1,2),f2],[(2,3),f3],[(3,10),f4]])
256
sage: f.intervals()
257
[(0, 1), (1, 2), (2, 3), (3, 10)]
258
"""
259
return self._intervals
260
261
def domain(self):
262
"""
263
Returns the domain of the function.
264
265
EXAMPLES::
266
267
sage: f1(x) = 1
268
sage: f2(x) = 1-x
269
sage: f3(x) = exp(x)
270
sage: f4(x) = sin(2*x)
271
sage: f = Piecewise([[(0,1),f1],[(1,2),f2],[(2,3),f3],[(3,10),f4]])
272
sage: f.domain()
273
(0, 10)
274
"""
275
endpoints = sum(self.intervals(), ())
276
return (min(endpoints), max(endpoints))
277
278
def functions(self):
279
"""
280
Returns the list of functions (the "pieces").
281
282
EXAMPLES::
283
284
sage: f1(x) = 1
285
sage: f2(x) = 1-x
286
sage: f3(x) = exp(x)
287
sage: f4(x) = sin(2*x)
288
sage: f = Piecewise([[(0,1),f1],[(1,2),f2],[(2,3),f3],[(3,10),f4]])
289
sage: f.functions()
290
[x |--> 1, x |--> -x + 1, x |--> e^x, x |--> sin(2*x)]
291
"""
292
return self._functions
293
294
def extend_by_zero_to(self,xmin=-1000,xmax=1000):
295
"""
296
This function simply returns the piecewise defined function which
297
is extended by 0 so it is defined on all of (xmin,xmax). This is
298
needed to add two piecewise functions in a reasonable way.
299
300
EXAMPLES::
301
302
sage: f1(x) = 1
303
sage: f2(x) = 1 - x
304
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
305
sage: f.extend_by_zero_to(-1, 3)
306
Piecewise defined function with 4 parts, [[(-1, 0), 0], [(0, 1), x |--> 1], [(1, 2), x |--> -x + 1], [(2, 3), 0]]
307
"""
308
zero = QQ['x'](0)
309
list_of_pairs = self.list()
310
a, b = self.domain()
311
if xmin < a:
312
list_of_pairs = [[(xmin, a), zero]] + list_of_pairs
313
if xmax > b:
314
list_of_pairs = list_of_pairs + [[(b, xmax), zero]]
315
return Piecewise(list_of_pairs)
316
317
def unextend(self):
318
"""
319
This removes any parts in the front or back of the function which
320
is zero (the inverse to extend_by_zero_to).
321
322
EXAMPLES::
323
324
sage: R.<x> = QQ[]
325
sage: f = Piecewise([[(-3,-1),1+2+x],[(-1,1),1-x^2]])
326
sage: e = f.extend_by_zero_to(-10,10); e
327
Piecewise defined function with 4 parts, [[(-10, -3), 0], [(-3, -1), x + 3], [(-1, 1), -x^2 + 1], [(1, 10), 0]]
328
sage: d = e.unextend(); d
329
Piecewise defined function with 2 parts, [[(-3, -1), x + 3], [(-1, 1), -x^2 + 1]]
330
sage: d==f
331
True
332
"""
333
list_of_pairs = self.list()
334
funcs = self.functions()
335
if funcs[0] == 0:
336
list_of_pairs = list_of_pairs[1:]
337
if funcs[-1] == 0:
338
list_of_pairs = list_of_pairs[:-1]
339
return Piecewise(list_of_pairs)
340
341
def _riemann_sum_helper(self, N, func, initial=0):
342
"""
343
A helper function for computing Riemann sums.
344
345
INPUT:
346
347
348
- ``N`` - the number of subdivisions
349
350
- ``func`` - a function to apply to the endpoints of
351
each subdivision
352
353
- ``initial`` - the starting value
354
355
356
EXAMPLES::
357
358
sage: f1(x) = x^2 ## example 1
359
sage: f2(x) = 5-x^2
360
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
361
sage: f._riemann_sum_helper(6, lambda x0, x1: (x1-x0)*f(x1))
362
19/6
363
"""
364
a,b = self.domain()
365
rsum = initial
366
h = (b-a)/N
367
for i in range(N):
368
x0 = a+i*h
369
x1 = a+(i+1)*h
370
rsum += func(x0, x1)
371
return rsum
372
373
def riemann_sum_integral_approximation(self,N,mode=None):
374
"""
375
Returns the piecewise line function defined by the Riemann sums in
376
numerical integration based on a subdivision into N subintervals.
377
378
Set mode="midpoint" for the height of the rectangles to be
379
determined by the midpoint of the subinterval; set mode="right" for
380
the height of the rectangles to be determined by the right-hand
381
endpoint of the subinterval; the default is mode="left" (the height
382
of the rectangles to be determined by the left-hand endpoint of
383
the subinterval).
384
385
EXAMPLES::
386
387
sage: f1(x) = x^2 ## example 1
388
sage: f2(x) = 5-x^2
389
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
390
sage: f.riemann_sum_integral_approximation(6)
391
17/6
392
sage: f.riemann_sum_integral_approximation(6,mode="right")
393
19/6
394
sage: f.riemann_sum_integral_approximation(6,mode="midpoint")
395
3
396
sage: f.integral(definite=True)
397
3
398
"""
399
if mode is None:
400
return self._riemann_sum_helper(N, lambda x0, x1: (x1-x0)*self(x0))
401
elif mode == "right":
402
return self._riemann_sum_helper(N, lambda x0, x1: (x1-x0)*self(x1))
403
elif mode == "midpoint":
404
return self._riemann_sum_helper(N, lambda x0, x1: (x1-x0)*self((x0+x1)/2))
405
else:
406
raise ValueError, "invalid mode"
407
408
def riemann_sum(self,N,mode=None):
409
"""
410
Returns the piecewise line function defined by the Riemann sums in
411
numerical integration based on a subdivision into N subintervals.
412
Set mode="midpoint" for the height of the rectangles to be
413
determined by the midpoint of the subinterval; set mode="right" for
414
the height of the rectangles to be determined by the right-hand
415
endpoint of the subinterval; the default is mode="left" (the height
416
of the rectangles to be determined by the left-hand endpoint of
417
the subinterval).
418
419
EXAMPLES::
420
421
sage: f1(x) = x^2
422
sage: f2(x) = 5-x^2
423
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
424
sage: f.riemann_sum(6,mode="midpoint")
425
Piecewise defined function with 6 parts, [[(0, 1/3), 1/36], [(1/3, 2/3), 1/4], [(2/3, 1), 25/36], [(1, 4/3), 131/36], [(4/3, 5/3), 11/4], [(5/3, 2), 59/36]]
426
427
::
428
429
sage: f = Piecewise([[(-1,1),(1-x^2).function(x)]])
430
sage: rsf = f.riemann_sum(7)
431
sage: P = f.plot(rgbcolor=(0.7,0.1,0.5), plot_points=40)
432
sage: Q = rsf.plot(rgbcolor=(0.7,0.6,0.6), plot_points=40)
433
sage: L = add([line([[a,0],[a,f(x=a)]],rgbcolor=(0.7,0.6,0.6)) for (a,b),f in rsf.list()])
434
sage: P + Q + L
435
436
::
437
438
sage: f = Piecewise([[(-1,1),(1/2+x-x^3)]], x) ## example 3
439
sage: rsf = f.riemann_sum(8)
440
sage: P = f.plot(rgbcolor=(0.7,0.1,0.5), plot_points=40)
441
sage: Q = rsf.plot(rgbcolor=(0.7,0.6,0.6), plot_points=40)
442
sage: L = add([line([[a,0],[a,f(x=a)]],rgbcolor=(0.7,0.6,0.6)) for (a,b),f in rsf.list()])
443
sage: P + Q + L
444
"""
445
if mode is None:
446
rsum = self._riemann_sum_helper(N, lambda x0,x1: [[(x0,x1),SR(self(x0))]],
447
initial=[])
448
elif mode == "right":
449
rsum = self._riemann_sum_helper(N, lambda x0,x1: [[(x0,x1),SR(self(x1))]],
450
initial=[])
451
elif mode == "midpoint":
452
rsum = self._riemann_sum_helper(N, lambda x0,x1: [[(x0,x1),SR(self((x0+x1)/2))]],
453
initial=[])
454
else:
455
raise ValueError, "invalid mode"
456
return Piecewise(rsum)
457
458
def trapezoid(self,N):
459
"""
460
Returns the piecewise line function defined by the trapezoid rule
461
for numerical integration based on a subdivision into N
462
subintervals.
463
464
EXAMPLES::
465
466
sage: R.<x> = QQ[]
467
sage: f1 = x^2
468
sage: f2 = 5-x^2
469
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
470
sage: f.trapezoid(4)
471
Piecewise defined function with 4 parts, [[(0, 1/2), 1/2*x], [(1/2, 1), 9/2*x - 2], [(1, 3/2), 1/2*x + 2], [(3/2, 2), -7/2*x + 8]]
472
473
::
474
475
sage: R.<x> = QQ[]
476
sage: f = Piecewise([[(-1,1),1-x^2]])
477
sage: tf = f.trapezoid(4)
478
sage: P = f.plot(rgbcolor=(0.7,0.1,0.5), plot_points=40)
479
sage: Q = tf.plot(rgbcolor=(0.7,0.6,0.6), plot_points=40)
480
sage: L = add([line([[a,0],[a,f(a)]],rgbcolor=(0.7,0.6,0.6)) for (a,b),f in tf.list()])
481
sage: P+Q+L
482
483
::
484
485
sage: R.<x> = QQ[]
486
sage: f = Piecewise([[(-1,1),1/2+x-x^3]]) ## example 3
487
sage: tf = f.trapezoid(6)
488
sage: P = f.plot(rgbcolor=(0.7,0.1,0.5), plot_points=40)
489
sage: Q = tf.plot(rgbcolor=(0.7,0.6,0.6), plot_points=40)
490
sage: L = add([line([[a,0],[a,f(a)]],rgbcolor=(0.7,0.6,0.6)) for (a,b),f in tf.list()])
491
sage: P+Q+L
492
"""
493
x = QQ['x'].gen()
494
def f(x0, x1):
495
f0, f1 = self(x0), self(x1)
496
return [[(x0,x1),f0+(f1-f0)*(x1-x0)**(-1)*(x-x0)]]
497
rsum = self._riemann_sum_helper(N, f, initial=[])
498
return Piecewise(rsum)
499
500
def trapezoid_integral_approximation(self,N):
501
"""
502
Returns the approximation given by the trapezoid rule for numerical
503
integration based on a subdivision into N subintervals.
504
505
EXAMPLES::
506
507
sage: f1(x) = x^2 ## example 1
508
sage: f2(x) = 1-(1-x)^2
509
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
510
sage: P = f.plot(rgbcolor=(0.7,0.1,0.5), plot_points=40)
511
sage: tf = f.trapezoid(6)
512
sage: Q = tf.plot(rgbcolor=(0.7,0.6,0.6), plot_points=40)
513
sage: ta = f.trapezoid_integral_approximation(6)
514
sage: t = text('trapezoid approximation = %s'%ta, (1.5, 0.25))
515
sage: a = f.integral(definite=True)
516
sage: tt = text('area under curve = %s'%a, (1.5, -0.5))
517
sage: P + Q + t + tt
518
519
::
520
521
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]]) ## example 2
522
sage: tf = f.trapezoid(4)
523
sage: ta = f.trapezoid_integral_approximation(4)
524
sage: Q = tf.plot(rgbcolor=(0.7,0.6,0.6), plot_points=40)
525
sage: t = text('trapezoid approximation = %s'%ta, (1.5, 0.25))
526
sage: a = f.integral(definite=True)
527
sage: tt = text('area under curve = %s'%a, (1.5, -0.5))
528
sage: P+Q+t+tt
529
"""
530
def f(x0, x1):
531
f0, f1 = self(x0), self(x1)
532
return ((f1+f0)/2)*(x1-x0)
533
return self._riemann_sum_helper(N, f)
534
535
def critical_points(self):
536
"""
537
Return the critical points of this piecewise function.
538
539
.. warning::
540
541
Uses maxima, which prints the warning to use results with
542
caution. Only works for piecewise functions whose parts are
543
polynomials with real critical not occurring on the
544
interval endpoints.
545
546
EXAMPLES::
547
548
sage: R.<x> = QQ[]
549
sage: f1 = x^0
550
sage: f2 = 10*x - x^2
551
sage: f3 = 3*x^4 - 156*x^3 + 3036*x^2 - 26208*x
552
sage: f = Piecewise([[(0,3),f1],[(3,10),f2],[(10,20),f3]])
553
sage: expected = [5, 12, 13, 14]
554
sage: all(abs(e-a) < 0.001 for e,a in zip(expected, f.critical_points()))
555
True
556
"""
557
from sage.calculus.calculus import maxima
558
x = var('x')
559
crit_pts = []
560
for (a,b), f in self.list():
561
for root in maxima.allroots(SR(f).diff(x)==0):
562
root = float(root.rhs())
563
if a < root < b:
564
crit_pts.append(root)
565
return crit_pts
566
567
def base_ring(self):
568
"""
569
Returns the base ring of the function pieces. This
570
is useful when this class is extended.
571
572
EXAMPLES::
573
574
sage: f1(x) = 1
575
sage: f2(x) = 1-x
576
sage: f3(x) = x^2-5
577
sage: f = Piecewise([[(0,1),f1],[(1,2),f2],[(2,3),f3]])
578
sage: base_ring(f)
579
Symbolic Ring
580
581
::
582
583
sage: R.<x> = QQ[]
584
sage: f1 = x^0
585
sage: f2 = 10*x - x^2
586
sage: f3 = 3*x^4 - 156*x^3 + 3036*x^2 - 26208*x
587
sage: f = Piecewise([[(0,3),f1],[(3,10),f2],[(10,20),f3]])
588
sage: f.base_ring()
589
Rational Field
590
"""
591
return (self.functions()[0]).base_ring()
592
593
def end_points(self):
594
"""
595
Returns a list of all interval endpoints for this function.
596
597
EXAMPLES::
598
599
sage: f1(x) = 1
600
sage: f2(x) = 1-x
601
sage: f3(x) = x^2-5
602
sage: f = Piecewise([[(0,1),f1],[(1,2),f2],[(2,3),f3]])
603
sage: f.end_points()
604
[0, 1, 2, 3]
605
"""
606
intervals = self.intervals()
607
return [ intervals[0][0] ] + [b for a,b in intervals]
608
609
def __call__(self,x0):
610
"""
611
Evaluates self at x0. Returns the average value of the jump if x0
612
is an interior endpoint of one of the intervals of self and the
613
usual value otherwise.
614
615
EXAMPLES::
616
617
sage: f1(x) = 1
618
sage: f2(x) = 1-x
619
sage: f3(x) = exp(x)
620
sage: f4(x) = sin(2*x)
621
sage: f = Piecewise([[(0,1),f1],[(1,2),f2],[(2,3),f3],[(3,10),f4]])
622
sage: f(0.5)
623
1
624
sage: f(5/2)
625
e^(5/2)
626
sage: f(5/2).n()
627
12.1824939607035
628
sage: f(1)
629
1/2
630
"""
631
#x0 = QQ(x0) ## does not allow for evaluation at pi
632
n = self.length()
633
endpts = self.end_points()
634
for i in range(1,n):
635
if x0 == endpts[i]:
636
return (self.functions()[i-1](x0) + self.functions()[i](x0))/2
637
if x0 == endpts[0]:
638
return self.functions()[0](x0)
639
if x0 == endpts[n]:
640
return self.functions()[n-1](x0)
641
for i in range(n):
642
if endpts[i] < x0 < endpts[i+1]:
643
return self.functions()[i](x0)
644
raise ValueError,"Value not defined outside of domain."
645
646
def which_function(self,x0):
647
"""
648
Returns the function piece used to evaluate self at x0.
649
650
EXAMPLES::
651
652
sage: f1(z) = z
653
sage: f2(x) = 1-x
654
sage: f3(y) = exp(y)
655
sage: f4(t) = sin(2*t)
656
sage: f = Piecewise([[(0,1),f1],[(1,2),f2],[(2,3),f3],[(3,10),f4]])
657
sage: f.which_function(3/2)
658
x |--> -x + 1
659
"""
660
for (a,b), f in self.list():
661
if a <= x0 <= b:
662
return f
663
raise ValueError,"Function not defined outside of domain."
664
665
def default_variable(self):
666
r"""
667
Return the default variable. The default variable is defined as the
668
first variable in the first piece uses a variable. If no pieces have
669
a variable (each piece is a constant value), `x` is returned.
670
671
The result is cached.
672
673
AUTHOR: Paul Butler
674
675
EXAMPLES::
676
677
sage: f1(x) = 1
678
sage: f2(x) = 5*x
679
sage: p = Piecewise([[(0,1),f1],[(1,4),f2]])
680
sage: p.default_variable()
681
x
682
683
sage: f1 = 3*var('y')
684
sage: p = Piecewise([[(0,1),4],[(1,4),f1]])
685
sage: p.default_variable()
686
y
687
688
"""
689
try:
690
return self.__default_variable
691
except AttributeError:
692
pass
693
for _, fun in self._list:
694
fun = SR(fun)
695
if fun.variables():
696
v = fun.variables()[0]
697
self.__default_variable = v
698
return v
699
# default to x
700
v = var('x')
701
self.__default_value = v
702
return v
703
704
def integral(self, x=None, a=None, b=None, definite=False):
705
r"""
706
By default, returns the indefinite integral of the function.
707
If definite=True is given, returns the definite integral.
708
709
AUTHOR:
710
711
- Paul Butler
712
713
EXAMPLES::
714
715
sage: f1(x) = 1-x
716
sage: f = Piecewise([[(0,1),1],[(1,2),f1]])
717
sage: f.integral(definite=True)
718
1/2
719
720
::
721
722
sage: f1(x) = -1
723
sage: f2(x) = 2
724
sage: f = Piecewise([[(0,pi/2),f1],[(pi/2,pi),f2]])
725
sage: f.integral(definite=True)
726
1/2*pi
727
728
sage: f1(x) = 2
729
sage: f2(x) = 3 - x
730
sage: f = Piecewise([[(-2, 0), f1], [(0, 3), f2]])
731
sage: f.integral()
732
Piecewise defined function with 2 parts, [[(-2, 0), x |--> 2*x + 4], [(0, 3), x |--> -1/2*x^2 + 3*x + 4]]
733
734
sage: f1(y) = -1
735
sage: f2(y) = y + 3
736
sage: f3(y) = -y - 1
737
sage: f4(y) = y^2 - 1
738
sage: f5(y) = 3
739
sage: f = Piecewise([[(-4,-3),f1],[(-3,-2),f2],[(-2,0),f3],[(0,2),f4],[(2,3),f5]])
740
sage: F = f.integral(y)
741
sage: F
742
Piecewise defined function with 5 parts, [[(-4, -3), y |--> -y - 4], [(-3, -2), y |--> 1/2*y^2 + 3*y + 7/2], [(-2, 0), y |--> -1/2*y^2 - y - 1/2], [(0, 2), y |--> 1/3*y^3 - y - 1/2], [(2, 3), y |--> 3*y - 35/6]]
743
744
Ensure results are consistent with FTC::
745
746
sage: F(-3) - F(-4)
747
-1
748
sage: F(-1) - F(-3)
749
1
750
sage: F(2) - F(0)
751
2/3
752
sage: f.integral(y, 0, 2)
753
2/3
754
sage: F(3) - F(-4)
755
19/6
756
sage: f.integral(y, -4, 3)
757
19/6
758
sage: f.integral(definite=True)
759
19/6
760
761
::
762
763
sage: f1(y) = (y+3)^2
764
sage: f2(y) = y+3
765
sage: f3(y) = 3
766
sage: f = Piecewise([[(-infinity, -3), f1], [(-3, 0), f2], [(0, infinity), f3]])
767
sage: f.integral()
768
Piecewise defined function with 3 parts, [[(-Infinity, -3), y |--> 1/3*y^3 + 3*y^2 + 9*y + 9], [(-3, 0), y |--> 1/2*y^2 + 3*y + 9/2], [(0, +Infinity), y |--> 3*y + 9/2]]
769
770
::
771
772
sage: f1(x) = e^(-abs(x))
773
sage: f = Piecewise([[(-infinity, infinity), f1]])
774
sage: f.integral(definite=True)
775
2
776
sage: f.integral()
777
Piecewise defined function with 1 parts, [[(-Infinity, +Infinity), x |--> -1/2*((sgn(x) - 1)*e^(2*x) - 2*e^x*sgn(x) + sgn(x) + 1)*e^(-x) - 1]]
778
779
::
780
781
sage: f = Piecewise([((0, 5), cos(x))])
782
sage: f.integral()
783
Piecewise defined function with 1 parts, [[(0, 5), x |--> sin(x)]]
784
785
786
TESTS:
787
788
Verify that piecewise integrals of zero work (trac #10841)::
789
790
sage: f0(x) = 0
791
sage: f = Piecewise([[(0,1),f0]])
792
sage: f.integral(x,0,1)
793
0
794
sage: f = Piecewise([[(0,1), 0]])
795
sage: f.integral(x,0,1)
796
0
797
sage: f = Piecewise([[(0,1), SR(0)]])
798
sage: f.integral(x,0,1)
799
0
800
801
"""
802
if a != None and b != None:
803
F = self.integral(x)
804
return F(b) - F(a)
805
806
if a != None or b != None:
807
raise TypeError, 'only one endpoint given'
808
809
area = 0 # cumulative definite integral of parts to the left of the current interval
810
integrand_pieces = self.list()
811
integrand_pieces.sort()
812
new_pieces = []
813
814
if x == None:
815
x = self.default_variable()
816
817
# The integral is computed by iterating over the pieces in order.
818
# The definite integral for each piece is calculated and accumulated in `area`.
819
# Thus at any time, `area` represents the definite integral of all the pieces
820
# encountered so far. The indefinite integral of each piece is also calculated,
821
# and the `area` before each piece is added to the piece.
822
#
823
# If a definite integral is requested, `area` is returned.
824
# Otherwise, a piecewise function is constructed from the indefinite integrals
825
# and returned.
826
#
827
# An exception is made if integral is called on a piecewise function
828
# that starts at -infinity. In this case, we do not try to calculate the
829
# definite integral of the first piece, and the value of `area` remains 0
830
# after the first piece.
831
832
for (start, end), fun in integrand_pieces:
833
fun = SR(fun)
834
if start == -infinity and not definite:
835
fun_integrated = fun.integral(x, end, x)
836
else:
837
try:
838
assume(start < x)
839
except ValueError: # Assumption is redundant
840
pass
841
fun_integrated = fun.integral(x, start, x) + area
842
forget(start < x)
843
if definite or end != infinity:
844
area += fun.integral(x, start, end)
845
new_pieces.append([(start, end), SR(fun_integrated).function(x)])
846
847
if definite:
848
return SR(area)
849
else:
850
return Piecewise(new_pieces)
851
852
def convolution(self, other):
853
"""
854
Returns the convolution function,
855
`f*g(t)=\int_{-\infty}^\infty f(u)g(t-u)du`, for compactly
856
supported `f,g`.
857
858
EXAMPLES::
859
860
sage: x = PolynomialRing(QQ,'x').gen()
861
sage: f = Piecewise([[(0,1),1*x^0]]) ## example 0
862
sage: g = f.convolution(f)
863
sage: h = f.convolution(g)
864
sage: P = f.plot(); Q = g.plot(rgbcolor=(1,1,0)); R = h.plot(rgbcolor=(0,1,1));
865
sage: # Type show(P+Q+R) to view
866
sage: f = Piecewise([[(0,1),1*x^0],[(1,2),2*x^0],[(2,3),1*x^0]]) ## example 1
867
sage: g = f.convolution(f)
868
sage: h = f.convolution(g)
869
sage: P = f.plot(); Q = g.plot(rgbcolor=(1,1,0)); R = h.plot(rgbcolor=(0,1,1));
870
sage: # Type show(P+Q+R) to view
871
sage: f = Piecewise([[(-1,1),1]]) ## example 2
872
sage: g = Piecewise([[(0,3),x]])
873
sage: f.convolution(g)
874
Piecewise defined function with 3 parts, [[(-1, 1), 0], [(1, 2), -3/2*x], [(2, 4), -3/2*x]]
875
sage: g = Piecewise([[(0,3),1*x^0],[(3,4),2*x^0]])
876
sage: f.convolution(g)
877
Piecewise defined function with 5 parts, [[(-1, 1), x + 1], [(1, 2), 3], [(2, 3), x], [(3, 4), -x + 8], [(4, 5), -2*x + 10]]
878
"""
879
f = self
880
g = other
881
M = min(min(f.end_points()),min(g.end_points()))
882
N = max(max(f.end_points()),max(g.end_points()))
883
R2 = PolynomialRing(QQ,2,names=["tt","uu"])
884
tt,uu = R2.gens()
885
conv = 0
886
f0 = f.functions()[0]
887
g0 = g.functions()[0]
888
R1 = f0.parent()
889
xx = R1.gen()
890
var = repr(xx)
891
if len(f.intervals())==1 and len(g.intervals())==1:
892
f = f.unextend()
893
g = g.unextend()
894
a1 = f.intervals()[0][0]
895
a2 = f.intervals()[0][1]
896
b1 = g.intervals()[0][0]
897
b2 = g.intervals()[0][1]
898
i1 = repr(f0).replace(var,repr(uu))
899
i2 = repr(g0).replace(var,"("+repr(tt-uu)+")")
900
cmd1 = "integrate((%s)*(%s),%s,%s,%s)"%(i1,i2, uu, a1, tt-b1) ## if a1+b1 < tt < a2+b1
901
cmd2 = "integrate((%s)*(%s),%s,%s,%s)"%(i1,i2, uu, tt-b2, tt-b1) ## if a1+b2 < tt < a2+b1
902
cmd3 = "integrate((%s)*(%s),%s,%s,%s)"%(i1,i2, uu, tt-b2, a2) ## if a1+b2 < tt < a2+b2
903
cmd4 = "integrate((%s)*(%s),%s,%s,%s)"%(i1,i2, uu, a1, a2) ## if a2+b1 < tt < a1+b2
904
conv1 = maxima.eval(cmd1)
905
conv2 = maxima.eval(cmd2)
906
conv3 = maxima.eval(cmd3)
907
conv4 = maxima.eval(cmd4)
908
# this is a very, very, very ugly hack
909
x = PolynomialRing(QQ,'x').gen()
910
fg1 = sage_eval(conv1.replace("tt",var), {'x':x}) ## should be = R2(conv1)
911
fg2 = sage_eval(conv2.replace("tt",var), {'x':x}) ## should be = R2(conv2)
912
fg3 = sage_eval(conv3.replace("tt",var), {'x':x}) ## should be = R2(conv3)
913
fg4 = sage_eval(conv4.replace("tt",var), {'x':x}) ## should be = R2(conv4)
914
if a1-b1<a2-b2:
915
if a2+b1!=a1+b2:
916
h = Piecewise([[(a1+b1,a1+b2),fg1],[(a1+b2,a2+b1),fg4],[(a2+b1,a2+b2),fg3]])
917
else:
918
h = Piecewise([[(a1+b1,a1+b2),fg1],[(a1+b2,a2+b2),fg3]])
919
else:
920
if a1+b2!=a2+b1:
921
h = Piecewise([[(a1+b1,a2+b1),fg1],[(a2+b1,a1+b2),fg2],[(a1+b2,a2+b2),fg3]])
922
else:
923
h = Piecewise([[(a1+b1,a2+b1),fg1],[(a2+b1,a2+b2),fg3]])
924
return h
925
926
if len(f.intervals())>1 or len(g.intervals())>1:
927
z = Piecewise([[(-3*abs(N-M),3*abs(N-M)),0*xx**0]])
928
ff = f.functions()
929
gg = g.functions()
930
intvlsf = f.intervals()
931
intvlsg = g.intervals()
932
for i in range(len(ff)):
933
for j in range(len(gg)):
934
f0 = Piecewise([[intvlsf[i],ff[i]]])
935
g0 = Piecewise([[intvlsg[j],gg[j]]])
936
h = g0.convolution(f0)
937
z = z + h
938
return z.unextend()
939
940
def derivative(self):
941
r"""
942
Returns the derivative (as computed by maxima)
943
Piecewise(I,`(d/dx)(self|_I)`), as I runs over the
944
intervals belonging to self. self must be piecewise polynomial.
945
946
EXAMPLES::
947
948
sage: f1(x) = 1
949
sage: f2(x) = 1-x
950
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
951
sage: f.derivative()
952
Piecewise defined function with 2 parts, [[(0, 1), x |--> 0], [(1, 2), x |--> -1]]
953
sage: f1(x) = -1
954
sage: f2(x) = 2
955
sage: f = Piecewise([[(0,pi/2),f1],[(pi/2,pi),f2]])
956
sage: f.derivative()
957
Piecewise defined function with 2 parts, [[(0, 1/2*pi), x |--> 0], [(1/2*pi, pi), x |--> 0]]
958
959
::
960
961
sage: f = Piecewise([[(0,1), (x * 2)]], x)
962
sage: f.derivative()
963
Piecewise defined function with 1 parts, [[(0, 1), x |--> 2]]
964
"""
965
x = var('x')
966
dlist = [[(a, b), derivative(f(x), x).function(x)] for (a,b),f in self.list()]
967
return Piecewise(dlist)
968
969
def tangent_line(self, pt):
970
"""
971
Computes the linear function defining the tangent line of the
972
piecewise function self.
973
974
EXAMPLES::
975
976
sage: f1(x) = x^2
977
sage: f2(x) = 5-x^3+x
978
sage: f = Piecewise([[(0,1),f1],[(1,2),f2]])
979
sage: tf = f.tangent_line(0.9) ## tangent line at x=0.9
980
sage: P = f.plot(rgbcolor=(0.7,0.1,0.5), plot_points=40)
981
sage: Q = tf.plot(rgbcolor=(0.7,0.2,0.2), plot_points=40)
982
sage: P + Q
983
"""
984
pt = QQ(pt)
985
R = QQ['x']
986
x = R.gen()
987
der = self.derivative()
988
tanline = (x-pt)*der(pt)+self(pt)
989
dlist = [[(a, b), tanline] for (a,b),f in self.list()]
990
return Piecewise(dlist)
991
992
def plot(self, *args, **kwds):
993
"""
994
Returns the plot of self.
995
996
Keyword arguments are passed onto the plot command for each piece
997
of the function. E.g., the plot_points keyword affects each
998
segment of the plot.
999
1000
EXAMPLES::
1001
1002
sage: f1(x) = 1
1003
sage: f2(x) = 1-x
1004
sage: f3(x) = exp(x)
1005
sage: f4(x) = sin(2*x)
1006
sage: f = Piecewise([[(0,1),f1],[(1,2),f2],[(2,3),f3],[(3,10),f4]])
1007
sage: P = f.plot(rgbcolor=(0.7,0.1,0), plot_points=40)
1008
sage: P
1009
1010
Remember: to view this, type show(P) or P.save("path/myplot.png")
1011
and then open it in a graphics viewer such as GIMP.
1012
1013
TESTS:
1014
1015
We should not add each piece to the legend individually, since
1016
this creates duplicates (:trac:`12651`). This tests that only
1017
one of the graphics objects in the plot has a non-``None``
1018
``legend_label``::
1019
1020
sage: f1 = sin(x)
1021
sage: f2 = cos(x)
1022
sage: f = piecewise([[(-1,0), f1],[(0,1), f2]])
1023
sage: p = f.plot(legend_label='$f(x)$')
1024
sage: lines = [
1025
... line
1026
... for line in p._Graphics__objects
1027
... if line.options()['legend_label'] is not None ]
1028
sage: len(lines)
1029
1
1030
"""
1031
from sage.plot.all import plot, Graphics
1032
1033
g = Graphics()
1034
1035
for i, ((a,b), f) in enumerate(self.list()):
1036
# If it's the first piece, pass all arguments. Otherwise,
1037
# filter out 'legend_label' so that we don't add each
1038
# piece to the legend separately (trac #12651).
1039
if i != 0 and 'legend_label' in kwds:
1040
del kwds['legend_label']
1041
1042
g += plot(f, a, b, *args, **kwds)
1043
1044
return g
1045
1046
def fourier_series_cosine_coefficient(self,n,L):
1047
r"""
1048
Returns the n-th Fourier series coefficient of
1049
`\cos(n\pi x/L)`, `a_n`.
1050
1051
INPUT:
1052
1053
1054
- ``self`` - the function f(x), defined over -L x L
1055
1056
- ``n`` - an integer n=0
1057
1058
- ``L`` - (the period)/2
1059
1060
1061
OUTPUT:
1062
`a_n = \frac{1}{L}\int_{-L}^L f(x)\cos(n\pi x/L)dx`
1063
1064
EXAMPLES::
1065
1066
sage: f(x) = x^2
1067
sage: f = Piecewise([[(-1,1),f]])
1068
sage: f.fourier_series_cosine_coefficient(2,1)
1069
pi^(-2)
1070
sage: f(x) = x^2
1071
sage: f = Piecewise([[(-pi,pi),f]])
1072
sage: f.fourier_series_cosine_coefficient(2,pi)
1073
1
1074
sage: f1(x) = -1
1075
sage: f2(x) = 2
1076
sage: f = Piecewise([[(-pi,pi/2),f1],[(pi/2,pi),f2]])
1077
sage: f.fourier_series_cosine_coefficient(5,pi)
1078
-3/5/pi
1079
"""
1080
from sage.all import cos, pi
1081
x = var('x')
1082
result = sum([(f(x)*cos(pi*x*n/L)/L).integrate(x, a, b)
1083
for (a,b), f in self.list()])
1084
if is_Expression(result):
1085
return result.simplify_trig()
1086
return result
1087
1088
def fourier_series_sine_coefficient(self,n,L):
1089
r"""
1090
Returns the n-th Fourier series coefficient of
1091
`\sin(n\pi x/L)`, `b_n`.
1092
1093
INPUT:
1094
1095
1096
- ``self`` - the function f(x), defined over -L x L
1097
1098
- ``n`` - an integer n0
1099
1100
- ``L`` - (the period)/2
1101
1102
1103
OUTPUT:
1104
`b_n = \frac{1}{L}\int_{-L}^L f(x)\sin(n\pi x/L)dx`
1105
1106
EXAMPLES::
1107
1108
sage: f(x) = x^2
1109
sage: f = Piecewise([[(-1,1),f]])
1110
sage: f.fourier_series_sine_coefficient(2,1) # L=1, n=2
1111
0
1112
"""
1113
from sage.all import sin, pi
1114
x = var('x')
1115
result = sum([(f(x)*sin(pi*x*n/L)/L).integrate(x, a, b)
1116
for (a,b), f in self.list()])
1117
if is_Expression(result):
1118
return result.simplify_trig()
1119
return result
1120
1121
def _fourier_series_helper(self, N, L, scale_function):
1122
r"""
1123
A helper function for the construction of Fourier series. The
1124
argument scale_function is a function which takes in n,
1125
representing the `n^{th}` coefficient, and return an
1126
expression to scale the sine and cosine coefficients by.
1127
1128
EXAMPLES::
1129
1130
sage: f(x) = x^2
1131
sage: f = Piecewise([[(-1,1),f]])
1132
sage: f._fourier_series_helper(3, 1, lambda n: 1)
1133
-4*cos(pi*x)/pi^2 + cos(2*pi*x)/pi^2 + 1/3
1134
"""
1135
from sage.all import pi, sin, cos, srange
1136
x = var('x')
1137
a0 = self.fourier_series_cosine_coefficient(0,L)
1138
result = a0/2 + sum([(self.fourier_series_cosine_coefficient(n,L)*cos(n*pi*x/L) +
1139
self.fourier_series_sine_coefficient(n,L)*sin(n*pi*x/L))*
1140
scale_function(n)
1141
for n in srange(1,N)])
1142
return result.expand()
1143
1144
1145
def fourier_series_partial_sum(self,N,L):
1146
r"""
1147
Returns the partial sum
1148
1149
.. math::
1150
1151
f(x) \sim \frac{a_0}{2} + \sum_{n=1}^N [a_n\cos(\frac{n\pi x}{L}) + b_n\sin(\frac{n\pi x}{L})],
1152
1153
as a string.
1154
1155
EXAMPLE::
1156
1157
sage: f(x) = x^2
1158
sage: f = Piecewise([[(-1,1),f]])
1159
sage: f.fourier_series_partial_sum(3,1)
1160
-4*cos(pi*x)/pi^2 + cos(2*pi*x)/pi^2 + 1/3
1161
sage: f1(x) = -1
1162
sage: f2(x) = 2
1163
sage: f = Piecewise([[(-pi,pi/2),f1],[(pi/2,pi),f2]])
1164
sage: f.fourier_series_partial_sum(3,pi)
1165
-3*sin(2*x)/pi + 3*sin(x)/pi - 3*cos(x)/pi - 1/4
1166
"""
1167
return self._fourier_series_helper(N, L, lambda n: 1)
1168
1169
def fourier_series_partial_sum_cesaro(self,N,L):
1170
r"""
1171
Returns the Cesaro partial sum
1172
1173
.. math::
1174
1175
f(x) \sim \frac{a_0}{2} + \sum_{n=1}^N (1-n/N)*[a_n\cos(\frac{n\pi x}{L}) + b_n\sin(\frac{n\pi x}{L})],
1176
1177
1178
as a string. This is a "smoother" partial sum - the Gibbs
1179
phenomenon is mollified.
1180
1181
EXAMPLE::
1182
1183
sage: f(x) = x^2
1184
sage: f = Piecewise([[(-1,1),f]])
1185
sage: f.fourier_series_partial_sum_cesaro(3,1)
1186
-8/3*cos(pi*x)/pi^2 + 1/3*cos(2*pi*x)/pi^2 + 1/3
1187
sage: f1(x) = -1
1188
sage: f2(x) = 2
1189
sage: f = Piecewise([[(-pi,pi/2),f1],[(pi/2,pi),f2]])
1190
sage: f.fourier_series_partial_sum_cesaro(3,pi)
1191
-sin(2*x)/pi + 2*sin(x)/pi - 2*cos(x)/pi - 1/4
1192
"""
1193
return self._fourier_series_helper(N, L, lambda n: 1-n/N)
1194
1195
def fourier_series_partial_sum_hann(self,N,L):
1196
r"""
1197
Returns the Hann-filtered partial sum (named after von Hann, not
1198
Hamming)
1199
1200
.. math::
1201
1202
f(x) \sim \frac{a_0}{2} + \sum_{n=1}^N H_N(n)*[a_n\cos(\frac{n\pi x}{L}) + b_n\sin(\frac{n\pi x}{L})],
1203
1204
as a string, where `H_N(x) = (1+\cos(\pi x/N))/2`. This is
1205
a "smoother" partial sum - the Gibbs phenomenon is mollified.
1206
1207
EXAMPLE::
1208
1209
sage: f(x) = x^2
1210
sage: f = Piecewise([[(-1,1),f]])
1211
sage: f.fourier_series_partial_sum_hann(3,1)
1212
-3*cos(pi*x)/pi^2 + 1/4*cos(2*pi*x)/pi^2 + 1/3
1213
sage: f1(x) = -1
1214
sage: f2(x) = 2
1215
sage: f = Piecewise([[(-pi,pi/2),f1],[(pi/2,pi),f2]])
1216
sage: f.fourier_series_partial_sum_hann(3,pi)
1217
-3/4*sin(2*x)/pi + 9/4*sin(x)/pi - 9/4*cos(x)/pi - 1/4
1218
"""
1219
from sage.all import cos, pi
1220
return self._fourier_series_helper(N, L, lambda n: (1+cos(pi*n/N))/2)
1221
1222
def fourier_series_partial_sum_filtered(self,N,L,F):
1223
r"""
1224
Returns the "filtered" partial sum
1225
1226
.. math::
1227
1228
f(x) \sim \frac{a_0}{2} + \sum_{n=1}^N F_n*[a_n\cos(\frac{n\pi x}{L}) + b_n\sin(\frac{n\pi x}{L})],
1229
1230
as a string, where `F = [F_1,F_2, ..., F_{N}]` is a list
1231
of length `N` consisting of real numbers. This can be used
1232
to plot FS solutions to the heat and wave PDEs.
1233
1234
EXAMPLE::
1235
1236
sage: f(x) = x^2
1237
sage: f = Piecewise([[(-1,1),f]])
1238
sage: f.fourier_series_partial_sum_filtered(3,1,[1,1,1])
1239
-4*cos(pi*x)/pi^2 + cos(2*pi*x)/pi^2 + 1/3
1240
sage: f1(x) = -1
1241
sage: f2(x) = 2
1242
sage: f = Piecewise([[(-pi,pi/2),f1],[(pi/2,pi),f2]])
1243
sage: f.fourier_series_partial_sum_filtered(3,pi,[1,1,1])
1244
-3*sin(2*x)/pi + 3*sin(x)/pi - 3*cos(x)/pi - 1/4
1245
"""
1246
return self._fourier_series_helper(N, L, lambda n: F[n])
1247
1248
def plot_fourier_series_partial_sum(self,N,L,xmin,xmax, **kwds):
1249
r"""
1250
Plots the partial sum
1251
1252
.. math::
1253
1254
f(x) \sim \frac{a_0}{2} + sum_{n=1}^N [a_n\cos(\frac{n\pi x}{L}) + b_n\sin(\frac{n\pi x}{L})],
1255
1256
over xmin x xmin.
1257
1258
EXAMPLE::
1259
1260
sage: f1(x) = -2
1261
sage: f2(x) = 1
1262
sage: f3(x) = -1
1263
sage: f4(x) = 2
1264
sage: f = Piecewise([[(-pi,-pi/2),f1],[(-pi/2,0),f2],[(0,pi/2),f3],[(pi/2,pi),f4]])
1265
sage: P = f.plot_fourier_series_partial_sum(3,pi,-5,5) # long time
1266
sage: f1(x) = -1
1267
sage: f2(x) = 2
1268
sage: f = Piecewise([[(-pi,pi/2),f1],[(pi/2,pi),f2]])
1269
sage: P = f.plot_fourier_series_partial_sum(15,pi,-5,5) # long time
1270
1271
Remember, to view this type show(P) or P.save("path/myplot.png")
1272
and then open it in a graphics viewer such as GIMP.
1273
"""
1274
from sage.plot.all import plot
1275
return plot(self.fourier_series_partial_sum(N,L), xmin, xmax, **kwds)
1276
1277
def plot_fourier_series_partial_sum_cesaro(self,N,L,xmin,xmax, **kwds):
1278
r"""
1279
Plots the partial sum
1280
1281
.. math::
1282
1283
f(x) \sim \frac{a_0}{2} + \sum_{n=1}^N (1-n/N)*[a_n\cos(\frac{n\pi x}{L}) + b_n\sin(\frac{n\pi x}{L})],
1284
1285
1286
over xmin x xmin. This is a "smoother" partial sum - the Gibbs
1287
phenomenon is mollified.
1288
1289
EXAMPLE::
1290
1291
sage: f1(x) = -2
1292
sage: f2(x) = 1
1293
sage: f3(x) = -1
1294
sage: f4(x) = 2
1295
sage: f = Piecewise([[(-pi,-pi/2),f1],[(-pi/2,0),f2],[(0,pi/2),f3],[(pi/2,pi),f4]])
1296
sage: P = f.plot_fourier_series_partial_sum_cesaro(3,pi,-5,5) # long time
1297
sage: f1(x) = -1
1298
sage: f2(x) = 2
1299
sage: f = Piecewise([[(-pi,pi/2),f1],[(pi/2,pi),f2]])
1300
sage: P = f.plot_fourier_series_partial_sum_cesaro(15,pi,-5,5) # long time
1301
1302
Remember, to view this type show(P) or P.save("path/myplot.png")
1303
and then open it in a graphics viewer such as GIMP.
1304
"""
1305
from sage.plot.all import plot
1306
return plot(self.fourier_series_partial_sum_cesaro(N,L), xmin, xmax, **kwds)
1307
1308
def plot_fourier_series_partial_sum_hann(self,N,L,xmin,xmax, **kwds):
1309
r"""
1310
Plots the partial sum
1311
1312
.. math::
1313
1314
f(x) \sim \frac{a_0}{2} + \sum_{n=1}^N H_N(n)*[a_n\cos(\frac{n\pi x}{L}) + b_n\sin(\frac{n\pi x}{L})],
1315
1316
1317
over xmin x xmin, where H_N(x) = (0.5)+(0.5)\*cos(x\*pi/N) is the
1318
N-th Hann filter.
1319
1320
EXAMPLE::
1321
1322
sage: f1(x) = -2
1323
sage: f2(x) = 1
1324
sage: f3(x) = -1
1325
sage: f4(x) = 2
1326
sage: f = Piecewise([[(-pi,-pi/2),f1],[(-pi/2,0),f2],[(0,pi/2),f3],[(pi/2,pi),f4]])
1327
sage: P = f.plot_fourier_series_partial_sum_hann(3,pi,-5,5) # long time
1328
sage: f1(x) = -1
1329
sage: f2(x) = 2
1330
sage: f = Piecewise([[(-pi,pi/2),f1],[(pi/2,pi),f2]])
1331
sage: P = f.plot_fourier_series_partial_sum_hann(15,pi,-5,5) # long time
1332
1333
Remember, to view this type show(P) or P.save("path/myplot.png")
1334
and then open it in a graphics viewer such as GIMP.
1335
"""
1336
from sage.plot.all import plot
1337
return plot(self.fourier_series_partial_sum_hann(N,L), xmin, xmax, **kwds)
1338
1339
def plot_fourier_series_partial_sum_filtered(self,N,L,F,xmin,xmax, **kwds):
1340
r"""
1341
Plots the partial sum
1342
1343
.. math::
1344
1345
f(x) \sim \frac{a_0}{2} + \sum_{n=1}^N F_n*[a_n\cos(\frac{n\pi x}{L}) + b_n\sin(\frac{n\pi x}{L})],
1346
1347
1348
over xmin x xmin, where `F = [F_1,F_2, ..., F_{N}]` is a
1349
list of length `N` consisting of real numbers. This can be
1350
used to plot FS solutions to the heat and wave PDEs.
1351
1352
EXAMPLE::
1353
1354
sage: f1(x) = -2
1355
sage: f2(x) = 1
1356
sage: f3(x) = -1
1357
sage: f4(x) = 2
1358
sage: f = Piecewise([[(-pi,-pi/2),f1],[(-pi/2,0),f2],[(0,pi/2),f3],[(pi/2,pi),f4]])
1359
sage: P = f.plot_fourier_series_partial_sum_filtered(3,pi,[1]*3,-5,5) # long time
1360
sage: f1(x) = -1
1361
sage: f2(x) = 2
1362
sage: f = Piecewise([[(-pi,-pi/2),f1],[(-pi/2,0),f2],[(0,pi/2),f1],[(pi/2,pi),f2]])
1363
sage: P = f.plot_fourier_series_partial_sum_filtered(15,pi,[1]*15,-5,5) # long time
1364
1365
Remember, to view this type show(P) or P.save("path/myplot.png")
1366
and then open it in a graphics viewer such as GIMP.
1367
"""
1368
from sage.plot.all import plot
1369
return plot(self.fourier_series_partial_sum_filtered(N,L,F), xmin, xmax, **kwds)
1370
1371
def fourier_series_value(self,x,L):
1372
r"""
1373
Returns the value of the Fourier series coefficient of self at
1374
`x`,
1375
1376
1377
.. math::
1378
1379
f(x) \sim \frac{a_0}{2} + \sum_{n=1}^\infty [a_n\cos(\frac{n\pi x}{L}) + b_n\sin(\frac{n\pi x}{L})], \ \ \ -L<x<L.
1380
1381
1382
This method applies to piecewise non-polynomial functions as well.
1383
1384
INPUT:
1385
1386
1387
- ``self`` - the function f(x), defined over -L x L
1388
1389
- ``x`` - a real number
1390
1391
- ``L`` - (the period)/2
1392
1393
1394
OUTPUT: `(f^*(x+)+f^*(x-)/2`, where `f^*` denotes
1395
the function `f` extended to `\RR` with period
1396
`2L` (Dirichlet's Theorem for Fourier series).
1397
1398
EXAMPLES::
1399
1400
sage: f1(x) = 1
1401
sage: f2(x) = 1-x
1402
sage: f3(x) = exp(x)
1403
sage: f4(x) = sin(2*x)
1404
sage: f = Piecewise([[(-10,1),f1],[(1,2),f2],[(2,3),f3],[(3,10),f4]])
1405
sage: f.fourier_series_value(101,10)
1406
1/2
1407
sage: f.fourier_series_value(100,10)
1408
1
1409
sage: f.fourier_series_value(10,10)
1410
1/2*sin(20) + 1/2
1411
sage: f.fourier_series_value(20,10)
1412
1
1413
sage: f.fourier_series_value(30,10)
1414
1/2*sin(20) + 1/2
1415
sage: f1(x) = -1
1416
sage: f2(x) = 2
1417
sage: f = Piecewise([[(-pi,0),lambda x:0],[(0,pi/2),f1],[(pi/2,pi),f2]])
1418
sage: f.fourier_series_value(-1,pi)
1419
0
1420
sage: f.fourier_series_value(20,pi)
1421
-1
1422
sage: f.fourier_series_value(pi/2,pi)
1423
1/2
1424
"""
1425
xnew = x - int(RR(x/(2*L)))*2*L
1426
endpts = self.end_points()
1427
if xnew == endpts[0] or xnew == endpts[-1]:
1428
return (self.functions()[0](endpts[0]) + self.functions()[-1](endpts[-1]))/2
1429
else:
1430
return self(xnew)
1431
1432
def cosine_series_coefficient(self,n,L):
1433
r"""
1434
Returns the n-th cosine series coefficient of
1435
`\cos(n\pi x/L)`, `a_n`.
1436
1437
INPUT:
1438
1439
1440
- ``self`` - the function f(x), defined over 0 x L (no
1441
checking is done to insure this)
1442
1443
- ``n`` - an integer n=0
1444
1445
- ``L`` - (the period)/2
1446
1447
1448
OUTPUT:
1449
`a_n = \frac{2}{L}\int_{-L}^L f(x)\cos(n\pi x/L)dx` such
1450
that
1451
1452
.. math::
1453
1454
f(x) \sim \frac{a_0}{2} + \sum_{n=1}^\infty a_n\cos(\frac{n\pi x}{L}),\ \ 0<x<L.
1455
1456
1457
1458
EXAMPLES::
1459
1460
sage: f(x) = x
1461
sage: f = Piecewise([[(0,1),f]])
1462
sage: f.cosine_series_coefficient(2,1)
1463
0
1464
sage: f.cosine_series_coefficient(3,1)
1465
-4/9/pi^2
1466
sage: f1(x) = -1
1467
sage: f2(x) = 2
1468
sage: f = Piecewise([[(0,pi/2),f1],[(pi/2,pi),f2]])
1469
sage: f.cosine_series_coefficient(2,pi)
1470
0
1471
sage: f.cosine_series_coefficient(3,pi)
1472
2/pi
1473
sage: f.cosine_series_coefficient(111,pi)
1474
2/37/pi
1475
sage: f1 = lambda x: x*(pi-x)
1476
sage: f = Piecewise([[(0,pi),f1]])
1477
sage: f.cosine_series_coefficient(0,pi)
1478
1/3*pi^2
1479
1480
"""
1481
from sage.all import cos, pi
1482
x = var('x')
1483
result = sum([(2*f(x)*cos(pi*x*n/L)/L).integrate(x, a, b)
1484
for (a,b), f in self.list()])
1485
if is_Expression(result):
1486
return result.simplify_trig()
1487
return result
1488
1489
1490
def sine_series_coefficient(self,n,L):
1491
r"""
1492
Returns the n-th sine series coefficient of
1493
`\sin(n\pi x/L)`, `b_n`.
1494
1495
INPUT:
1496
1497
- ``self`` - the function f(x), defined over 0 x L (no
1498
checking is done to insure this)
1499
1500
- ``n`` - an integer n0
1501
1502
- ``L`` - (the period)/2
1503
1504
OUTPUT:
1505
1506
`b_n = \frac{2}{L}\int_{-L}^L f(x)\sin(n\pi x/L)dx` such
1507
that
1508
1509
.. math::
1510
1511
f(x) \sim \sum_{n=1}^\infty b_n\sin(\frac{n\pi x}{L}),\ \ 0<x<L.
1512
1513
EXAMPLES::
1514
1515
sage: f(x) = 1
1516
sage: f = Piecewise([[(0,1),f]])
1517
sage: f.sine_series_coefficient(2,1)
1518
0
1519
sage: f.sine_series_coefficient(3,1)
1520
4/3/pi
1521
"""
1522
from sage.all import sin, pi
1523
x = var('x')
1524
result = sum([(2*f(x)*sin(pi*x*n/L)/L).integrate(x, a, b)
1525
for (a,b), f in self.list()])
1526
if is_Expression(result):
1527
return result.simplify_trig()
1528
return result
1529
1530
def laplace(self, x='x', s='t'):
1531
r"""
1532
Returns the Laplace transform of self with respect to the variable
1533
var.
1534
1535
INPUT:
1536
1537
1538
- ``x`` - variable of self
1539
1540
- ``s`` - variable of Laplace transform.
1541
1542
1543
We assume that a piecewise function is 0 outside of its domain and
1544
that the left-most endpoint of the domain is 0.
1545
1546
EXAMPLES::
1547
1548
sage: x, s, w = var('x, s, w')
1549
sage: f = Piecewise([[(0,1),1],[(1,2), 1-x]])
1550
sage: f.laplace(x, s)
1551
(s + 1)*e^(-2*s)/s^2 - e^(-s)/s + 1/s - e^(-s)/s^2
1552
sage: f.laplace(x, w)
1553
(w + 1)*e^(-2*w)/w^2 - e^(-w)/w + 1/w - e^(-w)/w^2
1554
1555
::
1556
1557
sage: y, t = var('y, t')
1558
sage: f = Piecewise([[(1,2), 1-y]])
1559
sage: f.laplace(y, t)
1560
(t + 1)*e^(-2*t)/t^2 - e^(-t)/t^2
1561
1562
::
1563
1564
sage: s = var('s')
1565
sage: t = var('t')
1566
sage: f1(t) = -t
1567
sage: f2(t) = 2
1568
sage: f = Piecewise([[(0,1),f1],[(1,infinity),f2]])
1569
sage: f.laplace(t,s)
1570
(s + 1)*e^(-s)/s^2 + 2*e^(-s)/s - 1/s^2
1571
"""
1572
from sage.all import assume, exp, forget
1573
x = var(x)
1574
s = var(s)
1575
assume(s>0)
1576
result = sum([(SR(f)*exp(-s*x)).integral(x,a,b)
1577
for (a,b),f in self.list()])
1578
forget(s>0)
1579
return result
1580
1581
def _make_compatible(self, other):
1582
"""
1583
Returns self and other extended to be defined on the same domain as
1584
well as a refinement of their intervals. This is used for adding
1585
and multiplying piecewise functions.
1586
1587
EXAMPLES::
1588
1589
sage: R.<x> = QQ[]
1590
sage: f1 = Piecewise([[(0, 2), x]])
1591
sage: f2 = Piecewise([[(1, 3), x^2]])
1592
sage: f1._make_compatible(f2)
1593
(Piecewise defined function with 2 parts, [[(0, 2), x], [(2, 3), 0]],
1594
Piecewise defined function with 2 parts, [[(0, 1), 0], [(1, 3), x^2]],
1595
[(0, 1), (1, 2), (2, 3)])
1596
"""
1597
a1, b1 = self.domain()
1598
a2, b2 = other.domain()
1599
a = min(a1, a2)
1600
b = max(b1, b2)
1601
F = self.extend_by_zero_to(a,b)
1602
G = other.extend_by_zero_to(a,b)
1603
endpts = list(set(F.end_points()).union(set(G.end_points())))
1604
endpts.sort()
1605
return F, G, zip(endpts, endpts[1:])
1606
1607
def __add__(self,other):
1608
"""
1609
Returns the piecewise defined function which is the sum of self and
1610
other. Does not require both domains be the same.
1611
1612
EXAMPLES::
1613
1614
sage: x = PolynomialRing(QQ,'x').gen()
1615
sage: f1 = x^0
1616
sage: f2 = 1-x
1617
sage: f3 = 2*x
1618
sage: f4 = 10-x
1619
sage: f = Piecewise([[(0,1),f1],[(1,2),f2],[(2,3),f3],[(3,10),f4]])
1620
sage: g1 = x-2
1621
sage: g2 = x-5
1622
sage: g = Piecewise([[(0,5),g1],[(5,10),g2]])
1623
sage: h = f+g
1624
sage: h
1625
Piecewise defined function with 5 parts, [[(0, 1), x - 1], [(1, 2), -1], [(2, 3), 3*x - 2], [(3, 5), 8], [(5, 10), 5]]
1626
1627
Note that in this case the functions must be defined using
1628
polynomial expressions *not* using the lambda notation.
1629
"""
1630
F, G, intervals = self._make_compatible(other)
1631
fcn = []
1632
for a,b in intervals:
1633
fcn.append([(a,b), F.which_function(b)+G.which_function(b)])
1634
return Piecewise(fcn)
1635
1636
def __mul__(self,other):
1637
r"""
1638
Returns the piecewise defined function which is the product of one
1639
piecewise function (self) with another one (other).
1640
1641
EXAMPLES::
1642
1643
sage: x = PolynomialRing(QQ,'x').gen()
1644
sage: f1 = x^0
1645
sage: f2 = 1-x
1646
sage: f3 = 2*x
1647
sage: f4 = 10-x
1648
sage: f = Piecewise([[(0,1),f1],[(1,2),f2],[(2,3),f3],[(3,10),f4]])
1649
sage: g1 = x-2
1650
sage: g2 = x-5
1651
sage: g = Piecewise([[(0,5),g1],[(5,10),g2]])
1652
sage: h = f*g
1653
sage: h
1654
Piecewise defined function with 5 parts, [[(0, 1), x - 2], [(1, 2), -x^2 + 3*x - 2], [(2, 3), 2*x^2 - 4*x], [(3, 5), -x^2 + 12*x - 20], [(5, 10), -x^2 + 15*x - 50]]
1655
sage: g*(11/2)
1656
Piecewise defined function with 2 parts, [[(0, 5), 11/2*x - 11], [(5, 10), 11/2*x - 55/2]]
1657
1658
Note that in this method the functions must be defined using
1659
polynomial expressions *not* using the lambda notation.
1660
"""
1661
## needed for scalar multiplication
1662
if isinstance(other,Rational) or isinstance(other,Integer):
1663
return Piecewise([[(a,b), other*f] for (a,b),f in self.list()])
1664
else:
1665
F, G, intervals = self._make_compatible(other)
1666
fcn = []
1667
for a,b in intervals:
1668
fcn.append([(a,b),F.which_function(b)*G.which_function(b)])
1669
return Piecewise(fcn)
1670
1671
__rmul__ = __mul__
1672
1673
def __eq__(self,other):
1674
"""
1675
Implements Boolean == operator.
1676
1677
EXAMPLES::
1678
1679
sage: f1 = x^0
1680
sage: f2 = 1-x
1681
sage: f3 = 2*x
1682
sage: f4 = 10-x
1683
sage: f = Piecewise([[(0,1),f1],[(1,2),f2],[(2,3),f3],[(3,10),f4]])
1684
sage: g = Piecewise([[(0,1),1],[(1,2),f2],[(2,3),f3],[(3,10),f4]])
1685
sage: f==g
1686
True
1687
"""
1688
return self.list()==other.list()
1689
1690