Path: blob/master/sage/plot/plot3d/parametric_surface.pyx
4045 views
r"""1Parametric Surface23Graphics 3D object for triangulating surfaces, and a base class for many other4objects that can be represented by a 2D parametrization.56It takes great care to turn degenerate quadrilaterals into triangles and7to propagate identified points to all attached polygons. This is not8so much to save space as it is to assist the raytracers/other rendering9systems to better understand the surface (and especially calculate correct10surface normals).1112AUTHORS:1314- Robert Bradshaw (2007-08-26): initial version1516EXAMPLES::1718sage: from sage.plot.plot3d.parametric_surface import ParametricSurface, MobiusStrip19sage: def f(x,y): return x+y, sin(x)*sin(y), x*y20sage: P = ParametricSurface(f, (srange(0,10,0.1), srange(-5,5.0,0.1)))21sage: show(P)22sage: S = MobiusStrip(1,.2)23sage: S.is_enclosed()24False25sage: S.show()2627.. NOTE::2829One may override ``eval()`` or ``eval_c()`` in a subclass30rather than passing in a function for greater speed.31One also would want to override get_grid.3233TODO: actually remove unused points, fix the below code::3435S = ParametricSurface(f=(lambda (x,y):(x,y,0)), domain=(range(10),range(10)))3637"""3839#*****************************************************************************40# Copyright (C) 2007 Robert Bradshaw <[email protected]>41#42# Distributed under the terms of the GNU General Public License (GPL)43#44# This code is distributed in the hope that it will be useful,45# but WITHOUT ANY WARRANTY; without even the implied warranty of46# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU47# General Public License for more details.48#49# The full text of the GPL is available at:50#51# http://www.gnu.org/licenses/52#*****************************************************************************53545556include "../../ext/stdsage.pxi"57include "../../ext/interrupt.pxi"5859include "point_c.pxi"6061from math import cos, sin62from sage.rings.all import RDF6364from base import RenderParams65from sage.ext.fast_eval cimport FastDoubleFunc66from sage.ext.interpreters.wrapper_rdf cimport Wrapper_rdf67from sage.ext.fast_eval import fast_float686970cdef inline bint smash_edge(point_c* vs, face_c* f, int a, int b):71if point_c_eq(vs[f.vertices[a]], vs[f.vertices[b]]):72f.vertices[b] = f.vertices[a]73f.n = 374return 175else:76return 07778cdef class ParametricSurface(IndexFaceSet):79"""80Base class that initializes the ParametricSurface81graphics type. This sets options, the function to be plotted, and the82plotting array as attributes.8384INPUT:8586- ``f`` - (default: None) The defining function. Either a tuple of87three functions, or a single function which returns a tuple, taking88two python floats as input. To subclass, pass None for f and override89eval_c or eval instead.9091- ``domain`` - (default: None) A tuple of two lists, defining the92grid of `u,v` values. If None, this will be calculate automatically.9394EXAMPLES::9596sage: from sage.plot.plot3d.parametric_surface import ParametricSurface97sage: def f(x,y): return cos(x)*sin(y), sin(x)*sin(y), cos(y)+log(tan(y/2))+0.2*x98sage: S = ParametricSurface(f, (srange(0,12.4,0.1), srange(0.1,2,0.1)))99sage: show(S)100101sage: len(S.face_list())1022214103104The Hessenberg surface:105106::107108sage: def f(u,v):109... a = 1110... from math import cos, sin, sinh, cosh111... x = cos(a)*(cos(u)*sinh(v)-cos(3*u)*sinh(3*v)/3) + sin(a)*(112... sin(u)*cosh(v)-sin(3*u)*cosh(3*v)/3)113... y = cos(a)*(sin(u)*sinh(v)+sin(3*u)*sinh(3*v)/3) + sin(a)*(114... -cos(u)*cosh(v)-cos(3*u)*cosh(3*v)/3)115... z = cos(a)*cos(2*u)*cosh(2*v)+sin(a)*sin(2*u)*sinh(2*v)116... return (x,y,z)117sage: v = srange(float(0),float((3/2)*pi),float(0.1))118sage: S = ParametricSurface(f, (srange(float(0),float(pi),float(0.1)),119... srange(float(-1),float(1),float(0.1))), color="blue")120sage: show(S)121"""122123def __init__(self, f=None, domain=None, **kwds):124"""125Create the graphics primitive :class:`ParametricSurface`. See the126docstring of this class for full documentation.127128EXAMPLES::129130sage: from sage.plot.plot3d.parametric_surface import ParametricSurface131sage: def f(x,y): return x+y, sin(x)*sin(y), x*y132sage: S = ParametricSurface(f, (srange(0,12.4,0.1), srange(0.1,2,0.1)))133"""134if isinstance(f, list):135f = tuple(f)136self.f = f137self.render_grid = domain138IndexFaceSet.__init__(self, [], [], **kwds)139140def default_render_params(self):141"""142Returns an instance of RenderParams suitable for plotting this object.143144TEST::145146sage: from sage.plot.plot3d.parametric_surface import MobiusStrip147sage: type(MobiusStrip(3,3).default_render_params())148<class 'sage.plot.plot3d.base.RenderParams'>149"""150return RenderParams(ds=.075, crease_threshold=.35)151152def x3d_geometry(self):153r"""154Returns XML-like representation of the coordinates of all points155in a triangulation of the object along with an indexing of those156points.157158TESTS::159160sage: _ = var('x,y')161sage: P = plot3d(x^2-y^2, (x, -2, 2), (y, -2, 2))162sage: s = P.x3d_str()163sage: s[:100]164"<Shape>\n<IndexedFaceSet coordIndex='0,1,..."165"""166self.triangulate(self.default_render_params())167return IndexFaceSet.x3d_geometry(self)168169def tachyon_repr(self, render_params):170"""171Returns representation of the object suitable for plotting172using Tachyon ray tracer.173174TESTS::175176sage: _ = var('x,y')177sage: P = plot3d(x^2-y^2, (x, -2, 2), (y, -2, 2))178sage: s = P.tachyon_repr(P.default_render_params())179sage: s[:2]180['TRI V0 -2 -2 0 V1 -2 -1.89744 0.399737 V2 -1.89744 -1.89744 0', 'texture...']181"""182self.triangulate(render_params)183return IndexFaceSet.tachyon_repr(self, render_params)184185def obj_repr(self, render_params):186"""187Returns complete representation of object with name, texture, and188lists of vertices, faces, and back-faces.189190TESTS::191192sage: _ = var('x,y')193sage: P = plot3d(x^2-y^2, (x, -2, 2), (y, -2, 2))194sage: s = P.obj_repr(P.default_render_params())195sage: s[:2]+s[2][:3]+s[3][:3]196['g obj_1', 'usemtl texture...', 'v -2 -2 0', 'v -2 -1.89744 0.399737', 'v -2 -1.79487 0.778435', 'f 1 2 42 41', 'f 2 3 43 42', 'f 3 4 44 43']197"""198self.triangulate(render_params)199return IndexFaceSet.obj_repr(self, render_params)200201def jmol_repr(self, render_params):202r"""203Returns representation of the object suitable for plotting204using Jmol.205206TESTS::207208sage: _ = var('x,y')209sage: P = plot3d(x^2-y^2, (x, -2, 2), (y, -2, 2))210sage: s = P.jmol_repr(P.testing_render_params())211sage: s[:10]212['pmesh obj_1 "obj_1.pmesh"\ncolor pmesh [102,102,255]']213"""214self.triangulate(render_params)215return IndexFaceSet.jmol_repr(self, render_params)216217def json_repr(self, render_params):218"""219Returns representation of the object in JSON format as220a list with one element, which is a string of a dictionary221listing vertices and faces.222223TESTS::224225sage: _ = var('x,y')226sage: P = plot3d(x^2-y^2, (x, -2, 2), (y, -2, 2))227sage: s = P.json_repr(P.default_render_params())228sage: s[0][:100]229'{vertices:[{x:-2,y:-2,z:0},{x:-2,y:-1.89744,z:0.399737},{x:-2,y:-1.79487,z:0.778435},{x:-2,y:-1.6923'230"""231self.triangulate(render_params)232return IndexFaceSet.json_repr(self, render_params)233234def is_enclosed(self):235"""236Returns a boolean telling whether or not it is necessary to237render the back sides of the polygons (assuming, of course,238that they have the correct orientation).239240This is calculated in by verifying the opposite edges241of the rendered domain either line up or are pinched together.242243EXAMPLES::244245sage: from sage.plot.plot3d.shapes import Sphere246sage: Sphere(1).is_enclosed()247True248249sage: from sage.plot.plot3d.parametric_surface import MobiusStrip250sage: MobiusStrip(1,0.2).is_enclosed()251False252"""253if self.fcount == 0:254self.triangulate()255return self.enclosed256257def dual(self):258"""259Returns an ``IndexFaceSet`` which is the dual of the260:class:`ParametricSurface` object as a triangulated surface.261262EXAMPLES:263264As one might expect, this gives an icosahedron::265266sage: D = dodecahedron()267sage: D.dual()268269But any enclosed surface should work::270271sage: from sage.plot.plot3d.shapes import Torus272sage: T = Torus(1, .2)273sage: T.dual()274sage: T.is_enclosed()275True276277Surfaces which are not enclosed, though, should raise an exception::278279sage: from sage.plot.plot3d.parametric_surface import MobiusStrip280sage: M = MobiusStrip(3,1)281sage: M.is_enclosed()282False283sage: M.dual()284Traceback (most recent call last):285...286NotImplementedError: This is only implemented for enclosed surfaces287288"""289# This doesn't completely make sense...290if self.fcount == 0:291self.triangulate()292if not self.is_enclosed():293raise NotImplementedError, "This is only implemented for enclosed surfaces"294return IndexFaceSet.dual(self)295296def bounding_box(self):297"""298Returns the lower and upper corners of a 3D bounding box for self.299This is used for rendering and self should fit entirely within this300box.301302Specifically, the first point returned should have x, y, and z303coordinates should be the respective infimum over all points in self,304and the second point is the supremum.305306EXAMPLES::307308sage: from sage.plot.plot3d.parametric_surface import MobiusStrip309sage: M = MobiusStrip(7,3,2)310sage: M.bounding_box()311((-10.0, -7.53907349250478..., -2.9940801852848145), (10.0, 7.53907349250478..., 2.9940801852848145))312"""313# We must triangulate before computing the bounding box; otherwise314# we'll get an empty bounding box, as the bounding box is computed315# using the triangulation, and before triangulating the triangulation316# is empty.317self.triangulate()318return IndexFaceSet.bounding_box(self)319320def triangulate(self, render_params=None):321r"""322Calls self.eval_grid() for all `(u,v)` in323`\text{urange} \times \text{vrange}` to construct this surface.324325The most complicated part of this code is identifying shared326vertices and shrinking trivial edges. This is not done so much327to save memory, rather it is needed so normals of the triangles328can be calculated correctly.329330TESTS::331332sage: from sage.plot.plot3d.parametric_surface import ParametricSurface, MobiusStrip333sage: def f(x,y): return x+y, sin(x)*sin(y), x*y # indirect doctests334sage: P = ParametricSurface(f, (srange(0,10,0.1), srange(-5,5.0,0.1))) # indirect doctests335sage: P.show() # indirect doctests336sage: S = MobiusStrip(1,.2) # indirect doctests337sage: S.show() # indirect doctests338"""339cdef double u, v340if render_params is None:341render_params = self.default_render_params()342ds = render_params.ds343if render_params.transform is not None:344ds /= render_params.transform.max_scale()345urange, vrange = self.get_grid(ds)346urange = [float(u) for u in urange]347vrange = [float(v) for v in vrange]348if self.render_grid == (urange, vrange) and self.fcount != 0:349# Already triangulated at on this grid.350return351352cdef Py_ssize_t i, j353cdef Py_ssize_t n = len(urange) - 1354cdef Py_ssize_t m = len(vrange) - 1355cdef Py_ssize_t ix = 0356357sig_on()358try:359self.realloc((m+1)*(n+1), m*n, 4*m*n)360self.eval_grid(urange, vrange)361except: # TODO -- this would catch control-C,etc. -- FIX THIS TO CATCH WHAT IS RAISED!!!!362sig_off()363self.fcount = self.vcount = 0364self.render_grid = None365raise366367# face_c.vertices:368#369# 0 - 1370# | |371# 3 - 2372373cdef face_c *face, *last_face374375for i from 0 <= i < n:376for j from 0 <= j < m:377ix = i*m+j378face = &self._faces[ix]379face.n = 4380face.vertices = &self.face_indices[4*ix]381382# Connect to the i-1 row383if i == 0:384if j == 0:385face.vertices[0] = 0386else:387face.vertices[0] = self._faces[ix-1].vertices[1]388face.vertices[1] = j+1389smash_edge(self.vs, face, 0, 1)390else:391face.vertices[0] = self._faces[ix-m].vertices[3]392face.vertices[1] = self._faces[ix-m].vertices[2]393394# Connect to the j-1 col395if j == 0:396face.vertices[3] = (i+1)*(m+1)397smash_edge(self.vs, face, 0, 3)398else:399face.vertices[3] = self._faces[ix-1].vertices[2]400401# This is the newly-seen vertex, identify if it's a triangle402face.vertices[2] = (i+1)*(m+1)+j+1403smash_edge(self.vs, face, 1, 2) or smash_edge(self.vs, face, 3, 2)404405# Now we see if it wraps around or is otherwise enclosed406cdef bint enclosed = 1407408cdef face_c *first, *last409for j from 0 <= j < m:410first = &self._faces[j]411last = &self._faces[(n-1)*m+j]412if point_c_eq(self.vs[first.vertices[0]], self.vs[last.vertices[3]]):413last.vertices[3] = first.vertices[0]414elif first.vertices[0] != first.vertices[1] or last.vertices[3] != last.vertices[2]:415enclosed = 0416if point_c_eq(self.vs[first.vertices[1]], self.vs[last.vertices[2]]):417last.vertices[2] = first.vertices[1]418elif first.vertices[0] != first.vertices[1] or last.vertices[3] != last.vertices[2]:419enclosed = 0420421for i from 0 <= i < n:422first = &self._faces[i*m]423last = &self._faces[i*m + m-1]424if point_c_eq(self.vs[first.vertices[0]], self.vs[last.vertices[1]]):425last.vertices[1] = first.vertices[0]426elif first.vertices[0] != first.vertices[3] or last.vertices[1] != last.vertices[2]:427enclosed = 0428if point_c_eq(self.vs[first.vertices[3]], self.vs[last.vertices[2]]):429last.vertices[2] = first.vertices[3]430elif first.vertices[0] != first.vertices[3] or last.vertices[1] != last.vertices[2]:431enclosed = 0432433self.enclosed = enclosed434435# make sure we deleted the correct point from the triangles436for ix from 0 <= ix < n*m:437face = &self._faces[ix]438if face.n == 3:439if face.vertices[3] == face.vertices[2] or face.vertices[3] == face.vertices[0]:440pass441else:442if face.vertices[0] == face.vertices[1]:443face.vertices[1] = face.vertices[2]444# face.vertices[1] == face.vertices[2]445face.vertices[2] = face.vertices[3]446447sig_off()448449self.vcount = (n+1)*(m+1)450self.fcount = n*m451self.icount = 4*n*m452self._clean_point_list()453454self.render_grid = urange, vrange455456457def get_grid(self, ds):458"""459TEST::460461sage: from sage.plot.plot3d.parametric_surface import ParametricSurface462sage: def f(x,y): return x+y,x-y,x*y463sage: P = ParametricSurface(f)464sage: P.get_grid(.1)465Traceback (most recent call last):466...467NotImplementedError: You must override the get_grid method.468"""469if self.render_grid is None:470raise NotImplementedError, "You must override the get_grid method."471return self.render_grid472473cdef int eval_grid(self, urange, vrange) except -1:474r"""475This fills in the points ``self.vs`` for all476`u \in \text{urange}, v \in \text{vrange}`.477We assume enough memory has been allocated.478479We branch outside the loops for efficiency. The options for self.f are:480481None -- call self.eval_c() or self.eval()482(One of these is presumably overridden.)483tuple -- split into fx, fy, fz and call each separately484callable -- call f(u,v)485486In addition, branches are taken for efficient calling of FastDoubleFunc487(including whether to iterate over python or c doubles).488"""489cdef Py_ssize_t i, j, m, n490cdef double u, v491cdef double uv[2]492cdef point_c *res493cdef double* ulist494cdef double* vlist495cdef bint fast_x, fast_y, fast_z496497if self.f is None:498499m, n = len(urange), len(vrange)500ulist = to_double_array(urange)501vlist = to_double_array(vrange)502503for i from 0 <= i < m:504u = ulist[i]505for j from 0 <= j < n:506v = vlist[j]507self.eval_c(&self.vs[i*n+j], u, v)508509sage_free(ulist)510sage_free(vlist)511512elif PY_TYPE_CHECK(self.f, tuple):513514fx, fy, fz = self.f515fast_x = PY_TYPE_CHECK(fx, FastDoubleFunc) or PY_TYPE_CHECK(fx, Wrapper_rdf)516fast_y = PY_TYPE_CHECK(fy, FastDoubleFunc) or PY_TYPE_CHECK(fx, Wrapper_rdf)517fast_z = PY_TYPE_CHECK(fz, FastDoubleFunc) or PY_TYPE_CHECK(fx, Wrapper_rdf)518519if fast_x or fast_y or fast_z:520521m, n = len(urange), len(vrange)522ulist = to_double_array(urange)523vlist = to_double_array(vrange)524525if PY_TYPE_CHECK(fx, FastDoubleFunc):526for i from 0 <= i < m:527uv[0] = ulist[i]528for j from 0 <= j < n:529uv[1] = vlist[j]530self.vs[i*n+j].x = (<FastDoubleFunc>fx)._call_c(uv)531elif fast_x: # must be Wrapper_rdf532for i from 0 <= i < m:533uv[0] = ulist[i]534for j from 0 <= j < n:535uv[1] = vlist[j]536(<Wrapper_rdf>fx).call_c(uv, &self.vs[i*n+j].x)537538539if PY_TYPE_CHECK(fy, FastDoubleFunc):540for i from 0 <= i < m:541uv[0] = ulist[i]542for j from 0 <= j < n:543uv[1] = vlist[j]544self.vs[i*n+j].y = (<FastDoubleFunc>fy)._call_c(uv)545elif fast_y: # must be Wrapper_rdf546for i from 0 <= i < m:547uv[0] = ulist[i]548for j from 0 <= j < n:549uv[1] = vlist[j]550(<Wrapper_rdf>fy).call_c(uv, &self.vs[i*n+j].y)551552if PY_TYPE_CHECK(fz, FastDoubleFunc):553for i from 0 <= i < m:554uv[0] = ulist[i]555for j from 0 <= j < n:556uv[1] = vlist[j]557self.vs[i*n+j].z = (<FastDoubleFunc>fz)._call_c(uv)558elif fast_z: # must be Wrapper_rdf559for i from 0 <= i < m:560uv[0] = ulist[i]561for j from 0 <= j < n:562uv[1] = vlist[j]563(<Wrapper_rdf>fz).call_c(uv, &self.vs[i*n+j].z)564565566sage_free(ulist)567sage_free(vlist)568569if not (fast_x and fast_y and fast_z):570ix = 0571for uu in urange:572for vv in vrange:573res = &self.vs[ix]574if not fast_x:575res.x = fx(uu, vv)576if not fast_y:577res.y = fy(uu, vv)578if not fast_z:579res.z = fz(uu, vv)580ix += 1581582else:583ix = 0584for uu in urange:585for vv in vrange:586res = &self.vs[ix]587res.x, res.y, res.z = self.f(uu, vv)588ix += 1589590# One of the following two methods should be overridden in591# derived classes.592593cdef int eval_c(self, point_c *res, double u, double v) except -1:594# can't do a cpdef because of the point_c* argument595res.x, res.y, res.z = self.eval(u, v)596597def eval(self, double u, double v):598"""599TEST::600601sage: from sage.plot.plot3d.parametric_surface import ParametricSurface602sage: def f(x,y): return x+y,x-y,x*y603sage: P = ParametricSurface(f,(srange(0,1,0.1),srange(0,1,0.1)))604sage: P.eval(0,0)605Traceback (most recent call last):606...607NotImplementedError608"""609raise NotImplementedError610611612class MobiusStrip(ParametricSurface):613"""614Base class for the :class:`MobiusStrip` graphics type. This sets the the615basic parameters of the object.616617INPUT:618619- ``r`` - A number which can be coerced to a float, serving roughly620as the radius of the object.621622- ``width`` - A number which can be coerced to a float, which gives the623width of the object.624625- ``twists`` - (default: 1) An integer, giving the number of twists in the626object (where one twist is the 'traditional' Mobius strip).627628EXAMPLES::629630sage: from sage.plot.plot3d.parametric_surface import MobiusStrip631sage: M = MobiusStrip(3,3)632sage: M.show()633"""634635def __init__(self, r, width, twists=1, **kwds):636"""637Create the graphics primitive MobiusStrip. See the docstring of638this class for full documentation.639640EXAMPLES:641642::643644sage: from sage.plot.plot3d.parametric_surface import MobiusStrip645sage: M = MobiusStrip(3,3); M # Same width and radius, roughly646sage: N = MobiusStrip(7,3,2); N # two twists, lots of open area in the middle647sage: O = MobiusStrip(5,1,plot_points=200,color='red'); O # keywords get passed to plot3d648649"""650ParametricSurface.__init__(self, **kwds)651self.r = float(r)652self.width = float(width)653self.twists = int(twists)654655def get_grid(self, ds):656"""657Returns appropriate `u` and `v` ranges for this MobiusStrip instance.658This is intended for internal use in creating an actual plot.659660INPUT:661662- ``ds`` - A number, typically coming from a RenderParams object,663which helps determine the increment for the `v` range for the664MobiusStrip object.665666EXAMPLE::667668sage: from sage.plot.plot3d.parametric_surface import MobiusStrip669sage: N = MobiusStrip(7,3,2) # two twists670sage: N.get_grid(N.default_render_params().ds)671([-1, 1], [0.0, 0.125663706144, 0.251327412287, 0.376991118431, ... 0])672673"""674twoPi = RDF.pi()*2675# Previous code, which doesn't seem to use any of the parameters676# TODO: figure out how to use it properly.677# res = max(min(twoPi*(self.r+self.twists*self.width)/ds, 10), 6*self.twists, 50)678res = max(6*self.twists, 50)679return [-1,1],[twoPi*k/res for k in range(res)] + [0]680681def eval(self, u, v):682"""683Returns tuple for `x,y,z` coordinates for the given ``u`` and ``v``684for this MobiusStrip instance.685686EXAMPLE::687688sage: from sage.plot.plot3d.parametric_surface import MobiusStrip689sage: N = MobiusStrip(7,3,2) # two twists690sage: N.eval(-1,0)691(4.0, 0.0, -0.0)692"""693return ( (self.r + u*self.width*cos(self.twists*v/2)) * cos(v),694(self.r + u*self.width*cos(self.twists*v/2)) * sin(v),695u*self.width*sin(self.twists*v/2) )696697698cdef double* to_double_array(py_list) except NULL:699cdef double* c_list = <double *>sage_malloc(sizeof(double) * len(py_list))700if c_list == NULL:701raise MemoryError702cdef Py_ssize_t i = 0703cdef double a704for a in py_list:705c_list[i] = a706i += 1707return c_list708709710711