Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

262

263

264

265

266

267

268

269

270

271

272

273

274

275

276

277

278

279

280

281

282

283

284

285

286

287

288

289

290

291

292

293

294

295

296

297

298

299

300

301

302

303

304

305

306

307

308

309

310

311

312

313

314

""" 

Complex Interpolation 

  

AUTHORS: 

  

- Ethan Van Andel (2009): initial version 

  

Development supported by NSF award No. 0702939. 

""" 

  

#***************************************************************************** 

# Copyright (C) 2009 Ethan Van Andel <evlutte@gmail.com>, 

# 

# Distributed under the terms of the GNU General Public License (GPL) 

# 

# This code is distributed in the hope that it will be useful, 

# but WITHOUT ANY WARRANTY; without even the implied warranty of 

# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 

# General Public License for more details. 

# 

# The full text of the GPL is available at: 

# 

# http://www.gnu.org/licenses/ 

#***************************************************************************** 

from __future__ import absolute_import 

  

import numpy as np 

cimport numpy as np 

  

from math import pi 

cdef double TWOPI = 2*pi 

  

def polygon_spline(pts): 

""" 

Creates a polygon from a set of complex or `(x,y)` points. The polygon 

will be a parametric curve from 0 to 2*pi. The returned values will be 

complex, not `(x,y)`. 

  

INPUT: 

  

- ``pts`` -- A list or array of complex numbers of tuples of the form 

`(x,y)`. 

  

EXAMPLES: 

  

A simple square:: 

  

sage: pts = [(-1, -1), (1, -1), (1, 1), (-1, 1)] 

sage: ps = polygon_spline(pts) 

sage: fx = lambda x: ps.value(x).real 

sage: fy = lambda x: ps.value(x).imag 

sage: show(parametric_plot((fx, fy), (0, 2*pi))) 

sage: m = Riemann_Map([lambda x: ps.value(real(x))], [lambda x: ps.derivative(real(x))],0) 

sage: show(m.plot_colored() + m.plot_spiderweb()) 

  

Polygon approximation of an circle:: 

  

sage: pts = [e^(I*t / 25) for t in range(25)] 

sage: ps = polygon_spline(pts) 

sage: ps.derivative(2) 

(-0.0470303661...+0.1520363883...j) 

""" 

return PSpline(pts) 

  

cdef class PSpline: 

""" 

A ``CCSpline`` object contains a polygon interpolation of a figure 

in the complex plane. 

  

EXAMPLES: 

  

A simple ``square``:: 

  

sage: pts = [(-1, -1), (1, -1), (1, 1), (-1, 1)] 

sage: ps = polygon_spline(pts) 

sage: ps.value(0) 

(-1-1j) 

sage: ps.derivative(0) 

(1.27323954...+0j) 

""" 

cdef int N 

cdef np.ndarray pts 

  

def __init__(self, pts): 

""" 

TESTS:: 

  

sage: pts = [(-1, -1), (1, -1), (1, 1), (-1, 1)] 

sage: ps = polygon_spline(pts) 

""" 

if type(pts[0]) == type((0,0)): 

self.pts = np.array( 

[np.complex(i[0], i[1]) for i in pts], dtype=np.complex128) 

else: 

self.pts = np.array(pts, dtype=np.complex128) 

self.N = len(pts) 

  

def value(self, double t): 

""" 

Returns the derivative (speed and direction of the curve) of a 

given point from the parameter ``t``. 

  

INPUT: 

  

- ``t`` -- double, the parameter value for the parameterized curve, 

between 0 and 2*pi. 

  

OUTPUT: 

  

A complex number representing the point on the polygon corresponding 

to the input ``t``. 

  

EXAMPLES: 

  

:: 

  

sage: pts = [(-1, -1), (1, -1), (1, 1), (-1, 1)] 

sage: ps = polygon_spline(pts) 

sage: ps.value(.5) 

(-0.363380227632...-1j) 

sage: ps.value(0) - ps.value(2*pi) 

0j 

sage: ps.value(10) 

(0.26760455264...+1j) 

""" 

cdef double t1 = ((t / TWOPI*self.N) % self.N) 

pt1 = self.pts[int(t1)] 

pt2 = self.pts[(int(t1) + 1) % self.N] 

return pt1 + (pt2 - pt1) * (t1 - int(t1)) 

  

def derivative(self, double t): 

""" 

Returns the derivative (speed and direction of the curve) of a 

given point from the parameter ``t``. 

  

INPUT: 

  

- ``t`` -- double, the parameter value for the parameterized curve, 

between 0 and 2*pi. 

  

OUTPUT: 

  

A complex number representing the derivative at the point on the 

polygon corresponding to the input ``t``. 

  

EXAMPLES: 

  

:: 

  

sage: pts = [(-1, -1), (1, -1), (1, 1), (-1, 1)] 

sage: ps = polygon_spline(pts) 

sage: ps.derivative(1 / 3) 

(1.27323954473...+0j) 

sage: ps.derivative(0) - ps.derivative(2*pi) 

0j 

sage: ps.derivative(10) 

(-1.27323954473...+0j) 

""" 

cdef double t1 = ((t / TWOPI*self.N) % self.N) 

pt1 = self.pts[int(t1)] 

pt2 = self.pts[(int(t1) + 1) % self.N] 

return (pt2 - pt1) * self.N / TWOPI 

  

def complex_cubic_spline(pts): 

""" 

Creates a cubic spline interpolated figure from a set of complex or 

`(x,y)` points. The figure will be a parametric curve from 0 to 2*pi. 

The returned values will be complex, not `(x,y)`. 

  

INPUT: 

  

- ``pts`` A list or array of complex numbers, or tuples of the form 

`(x,y)`. 

  

EXAMPLES: 

  

A simple ``square``:: 

  

sage: pts = [(-1, -1), (1, -1), (1, 1), (-1, 1)] 

sage: cs = complex_cubic_spline(pts) 

sage: fx = lambda x: cs.value(x).real 

sage: fy = lambda x: cs.value(x).imag 

sage: show(parametric_plot((fx, fy), (0, 2*pi))) 

sage: m = Riemann_Map([lambda x: cs.value(real(x))], [lambda x: cs.derivative(real(x))], 0) 

sage: show(m.plot_colored() + m.plot_spiderweb()) 

  

Polygon approximation of a circle:: 

  

sage: pts = [e^(I*t / 25) for t in range(25)] 

sage: cs = complex_cubic_spline(pts) 

sage: cs.derivative(2) 

(-0.0497765406583...+0.151095006434...j) 

""" 

return CCSpline(pts) 

  

cdef class CCSpline: 

""" 

A ``CCSpline`` object contains a cubic interpolation of a figure 

in the complex plane. 

  

EXAMPLES: 

  

A simple ``square``:: 

  

sage: pts = [(-1, -1), (1, -1), (1, 1), (-1, 1)] 

sage: cs = complex_cubic_spline(pts) 

sage: cs.value(0) 

(-1-1j) 

sage: cs.derivative(0) 

(0.9549296...-0.9549296...j) 

""" 

cdef int N 

cdef np.ndarray avec,bvec,cvec,dvec 

  

#standard cubic interpolation method 

def __init__(self, pts): 

""" 

TESTS:: 

  

sage: pts = [(-1, -1), (1, -1), (1, 1), (-1, 1)] 

sage: cs = complex_cubic_spline(pts) 

""" 

if type(pts[0]) == type((0,0)): 

pts = np.array( 

[np.complex(pt[0], pt[1]) for pt in pts], dtype=np.complex128) 

cdef int N, i, k 

N = len(pts) 

yvec = np.zeros(N, dtype=np.complex128) 

for i in xrange(N): 

yvec[i] = 3 * (pts[(i - 1) % N] - 2*pts[i] + pts[(i + 1) % N]) 

bmat = np.zeros([N, N], dtype=np.complex128) 

for i in xrange(N): 

bmat[i, i] = 4 

bmat[(i - 1) % N, i] = 1 

bmat[(i + 1) % N, i] = 1 

bvec = (np.linalg.solve(bmat, yvec)) 

cvec = np.zeros(N, dtype=np.complex128) 

for i in xrange(N): 

cvec[i] = (pts[(i + 1) % N] - pts[i] - 1.0/3.0 * 

bvec[(i + 1) % N] - 2./3. * bvec[i]) 

dvec = np.array(pts, dtype=np.complex128) 

avec = np.zeros(N, dtype=np.complex128) 

for i in xrange(N): 

avec[i] = 1.0/3.0 * (bvec[(i + 1) % N] - bvec[i]) 

self.avec = avec 

self.bvec = bvec 

self.cvec = cvec 

self.dvec = dvec 

self.N = N 

  

def value(self, double t): 

""" 

Returns the location of a given point from the parameter ``t``. 

  

INPUT: 

  

- ``t`` -- double, the parameter value for the parameterized curve, 

between 0 and 2*pi. 

  

OUTPUT: 

  

A complex number representing the point on the figure corresponding 

to the input ``t``. 

  

EXAMPLES: 

  

:: 

  

sage: pts = [(-1, -1), (1, -1), (1, 1), (-1, 1)] 

sage: cs = complex_cubic_spline(pts) 

sage: cs.value(4 / 7) 

(-0.303961332787...-1.34716728183...j) 

sage: cs.value(0) - cs.value(2*pi) 

0j 

sage: cs.value(-2.73452) 

(0.934561222231...+0.881366116402...j) 

""" 

cdef double t1 = (t / TWOPI * self.N) % self.N 

cdef int j = int(t1) 

return (self.avec[j] * (t1 - j)**3 + self.bvec[j] * (t1 - j)**2 + 

self.cvec[j] * (t1 - j) + self.dvec[j]) 

  

def derivative(self, double t): 

""" 

Returns the derivative (speed and direction of the curve) of a 

given point from the parameter ``t``. 

  

INPUT: 

  

- ``t`` -- double, the parameter value for the parameterized curve, 

between 0 and 2*pi. 

  

OUTPUT: 

  

A complex number representing the derivative at the point on the 

figure corresponding to the input ``t``. 

  

EXAMPLES: 

  

:: 

  

sage: pts = [(-1, -1), (1, -1), (1, 1), (-1, 1)] 

sage: cs = complex_cubic_spline(pts) 

sage: cs.derivative(3 / 5) 

(1.40578892327...-0.225417136326...j) 

sage: cs.derivative(0) - cs.derivative(2 * pi) 

0j 

sage: cs.derivative(-6) 

(2.52047692949...-1.89392588310...j) 

""" 

cdef double t1 = (t / TWOPI * self.N) % self.N 

cdef int j = int(t1) 

return (3 * self.avec[j] * (t1 - j)**2 + 2 * self.bvec[j] * (t1 - j) + 

self.cvec[j]) / TWOPI * self.N