Path: blob/master/sage/plot/plot3d/implicit_surface.pyx
4036 views
r"""1Graphics 3D object for representing and triangulating isosurfaces.23AUTHORS:45- Robert Hanson (2007): initial Java version, in Jmol.6- Carl Witty (2009-01): first Cython version.7- Bill Cauchois (2009): improvements for inclusion into Sage.8"""910#*****************************************************************************11# Copyright (C) 2009 Carl Witty <[email protected]>12#13# Distributed under the terms of the GNU General Public License (GPL)14#15# This code is distributed in the hope that it will be useful,16# but WITHOUT ANY WARRANTY; without even the implied warranty of17# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU18# General Public License for more details.19#20# The full text of the GPL is available at:21#22# http://www.gnu.org/licenses/23#*****************************************************************************2425# Pieces of this file are strongly based on the marching cubes26# implementation in Jmol located at src/org/jmol/jvxl/calc/MarchingCubes.java.27# The data tables are in fact directly copied from there.2829#*****************************************************************************30# This copyright is inherited from MarchingCubes.java in the Jmol source code.31#32# * Copyright (C) 2007 Miguel, Bob, Jmol Development33# *34# * Contact: [email protected]35# *36# * This library is free software; you can redistribute it and/or37# * modify it under the terms of the GNU Lesser General Public38# * License as published by the Free Software Foundation; either39# * version 2.1 of the License, or (at your option) any later version.40# *41# * This library is distributed in the hope that it will be useful,42# * but WITHOUT ANY WARRANTY; without even the implied warranty of43# * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU44# * Lesser General License for more details.45# *46# * You should have received a copy of the GNU Lesser General Public47# * License along with this library; if not, write to the Free Software48# * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA.49#*****************************************************************************505152# There's a framework in here for computing multiple isosurfaces of a53# single function. Currently, it's only used for a single54# implicit_plot3d where contour=... is given a list, but it's ready to55# be extended. I think the best UI would be if animate(...) and56# show(...) had a prepass where they went through their rendering57# trees, found ImplicitSurface objects with the same function,58# bounding box, and plot_points (other parameters, such as contours,59# hole, jmol_color, vertex_color, would be allowed to be different),60# and arranged to have them all plotted simultaneously. These61# prepasses have not been written yet. Combining multiple62# ImplicitSurface plots would be particularly advantageous for animate(...),63# but for a big animation, watch out for memory usage.6465# If you have a reasonably simple surface (not a space-filling fractal),66# then if n is your resolution, we have n^3 evaluations of the main67# function, about n^2 evaluations of auxiliary functions (hole, gradient,68# vertex_color/jmol_color), and output of size about n^2.6970# With this in mind, we pay particular attention to optimizing the n^371# function evaluations. (But keep in mind that n may be as small as 2072# or so, so we shouldn't ignore the efficiency of the n^2 parts.)7374# Also, we make sure not to ever allocate O(n^3) memory; we do the75# computation a slice at a time to avoid this. (Jmol always allocates76# n^3 memory when it reads the JVXL file, but that might be on a different77# computer; Tachyon would only allocate memory proportional to the78# output size.)7980from cStringIO import StringIO8182cimport numpy as np83import numpy as np8485from sage.plot.plot3d.transform cimport point_c, face_c, Transformation86from sage.plot.plot3d.base cimport PrimitiveObject87from sage.plot.plot3d.base import RenderParams, default_texture88from sage.plot.plot3d.index_face_set cimport IndexFaceSet89from sage.rings.all import RDF90from sage.plot.misc import setup_for_eval_on_grid9192include '../../ext/cdefs.pxi'93include '../../ext/stdsage.pxi'94include '../../gsl/gsl.pxi'95include "../../ext/python_string.pxi"9697include "point_c.pxi"9899# The default value for plot_points (a keyword argument to implicit_plot3d()),100# assumed when plot_points is set to "automatic".101DEFAULT_PLOT_POINTS = 40102103cdef double nan = float(RDF('NaN'))104105cdef inline bint marching_has_edge(double a, double b, double contour, double *f, bint *has_nan):106# XXX Would be nicer to use isnan(), because it's inlined.107# Is it portable enough?108if gsl_isnan(a) or gsl_isnan(b):109has_nan[0] = True110return False111112has_nan[0] = False113114if (a >= contour) == (b >= contour):115return False116117f[0] = (contour - a) / (b - a)118return True119120# Returns 0 or 1121cdef inline int marching_is_inside(double v, double contour):122return gsl_isnan(v) or v < contour123124cdef void interpolate_point_c(point_c *result, double frac, point_c *inputs):125result[0].x = inputs[0].x + frac*(inputs[1].x - inputs[0].x)126result[0].y = inputs[0].y + frac*(inputs[1].y - inputs[0].y)127result[0].z = inputs[0].z + frac*(inputs[1].z - inputs[0].z)128129cdef class VertexInfo:130# The point in "integer space"131cdef point_c pt132133# The gradient at this point in "evaluation space"134cdef point_c gradient135136# (R,G,B), so not really a point at all137cdef point_c color138139# This point in "evaluation space"140cdef point_c eval_pt141142cdef void update_eval_pt(self, point_c *eval_min, point_c *eval_scale):143"""144Use eval_min and eval_scale to transform self.pt into evaluation space145and store the result in self.eval_pt.146"""147self.eval_pt.x = eval_min[0].x + eval_scale[0].x*self.pt.x148self.eval_pt.y = eval_min[0].y + eval_scale[0].y*self.pt.y149self.eval_pt.z = eval_min[0].z + eval_scale[0].z*self.pt.z150151def __repr__(self):152"""153TESTS::154155sage: from sage.plot.plot3d.implicit_surface import VertexInfo156sage: VertexInfo()157<0.0, 0.0, 0.0>158"""159return '<%s, %s, %s>' % (self.pt.x, self.pt.y, self.pt.z)160161cdef mk_VertexInfo(double x, double y, double z, point_c *eval_min, point_c *eval_scale):162cdef VertexInfo v163v = PY_NEW(VertexInfo)164v.pt.x = x165v.pt.y = y166v.pt.z = z167168v.update_eval_pt(eval_min, eval_scale)169170return v171172cdef class MarchingCubes:173r"""174Handles marching cube rendering.175176Protocol:1771781. Create the class.1792. Call process_slice once for each X slice, from self.nx > x >= 0.1803. Call finish(), which returns a list of strings.181182Note: Actually, only 4 slices ever exist; the caller will re-use old183storage.184"""185186cdef readonly object xrange187cdef readonly object yrange188cdef readonly object zrange189cdef readonly double contour190cdef readonly int nx191cdef readonly int ny192cdef readonly int nz193cdef readonly Transformation transform194cdef readonly object region195cdef readonly object gradient196cdef readonly bint smooth197cdef readonly object vertex_color198cdef readonly object results199200# We deal with three coordinate systems. We do most of our work201# in an integral coordinate system where x ranges from202# 0 <= x < self.nx, etc.; we do function evaluations where203# self.xrange[0] <= x <= self.xrange[1], etc.; and output204# is in a third coordinate system.205# (Note that in "integer space", function evaluations of the main206# function happen on integer coordinates; but function evaluations207# of the auxiliary functions will have one non-integer coordinate.)208209# These parameters convert from integer space to function-evaluation210# space: eval_coord = eval_min + int_coord * eval_scale211cdef point_c eval_min212cdef point_c eval_scale213214# The componentwise reciprocal of eval_scale; just used to change215# some divisions into multiplications216cdef point_c eval_scale_inv217218cdef point_c out_origin, out_plus_x, out_plus_y, out_plus_z219220def __init__(self, xrange, yrange, zrange, contour, plot_points,221transform=None, region=None, gradient=None, smooth=True, vertex_color=None):222"""223TESTS:224225Marching cubes is an abstract base class, you can't do anything with it::226227sage: from sage.plot.plot3d.implicit_surface import MarchingCubes228sage: cube_marcher = MarchingCubes((0, 1), (0, 1), (0, 1), 1, (10, 10, 10))229"""230231self.xrange = xrange232self.yrange = yrange233self.zrange = zrange234self.contour = contour235self.nx = plot_points[0]236self.ny = plot_points[1]237self.nz = plot_points[2]238self.transform = transform239self.region = region240self.gradient = gradient241self.smooth = smooth242self.vertex_color = vertex_color243244self.eval_min.x = xrange[0]245self.eval_scale.x = (xrange[1] - xrange[0]) / (self.nx - 1)246self.eval_min.y = yrange[0]247self.eval_scale.y = (yrange[1] - yrange[0]) / (self.ny - 1)248self.eval_min.z = zrange[0]249self.eval_scale.z = (zrange[1] - zrange[0]) / (self.nz - 1)250self.eval_scale_inv.x = 1/self.eval_scale.x251self.eval_scale_inv.y = 1/self.eval_scale.y252self.eval_scale_inv.z = 1/self.eval_scale.z253254cdef point_c zero_pt, origin, plus_x, plus_y, plus_z255zero_pt.x = 0; zero_pt.y = 0; zero_pt.z = 0256origin = self.eval_min257plus_x = zero_pt; plus_x.x = self.eval_scale.x258plus_y = zero_pt; plus_y.y = self.eval_scale.y259plus_z = zero_pt; plus_z.z = self.eval_scale.z260if self.transform is not None:261self.transform.transform_point_c(&self.out_origin, origin)262self.transform.transform_point_c(&self.out_plus_x, plus_x)263self.transform.transform_point_c(&self.out_plus_y, plus_y)264self.transform.transform_point_c(&self.out_plus_z, plus_z)265else:266self.out_origin = origin267self.out_plus_x = plus_x268self.out_plus_y = plus_y269self.out_plus_z = plus_z270271self.results = []272273def finish(self):274"""275Returns the results of the marching cubes algorithm as a list. The format276is specific to the subclass implementing this method.277278TESTS:279280By default, it returns an empty list::281282sage: from sage.plot.plot3d.implicit_surface import MarchingCubes283sage: cube_marcher = MarchingCubes((0, 1), (0, 1), (0, 1), 1, (10, 10, 10), None)284sage: cube_marcher.finish()285[]286"""287return self.results288289cdef class MarchingCubesTriangles(MarchingCubes):290"""291A subclass of MarchingCubes that returns its results as a list of triangles,292including their vertices and normals (if smooth=True).293"""294295cdef readonly np.ndarray x_vertices296cdef readonly np.ndarray y_vertices297cdef readonly np.ndarray z_vertices298299cdef readonly np.ndarray y_vertices_swapped300cdef readonly np.ndarray z_vertices_swapped301302cdef readonly slices303304def __init__(self, *args, **kwargs):305"""306TESTS::307308sage: from sage.plot.plot3d.implicit_surface import MarchingCubesTriangles309sage: cube_marcher = MarchingCubesTriangles((0, 1), (0, 1), (0, 1), 1, (10, 10, 10))310"""311312MarchingCubes.__init__(self, *args, **kwargs)313314self.x_vertices = np.empty((self.ny, self.nz), dtype=object)315self.y_vertices = np.empty((2, self.ny-1, self.nz), dtype=object)316self.z_vertices = np.empty((2, self.ny, self.nz-1), dtype=object)317318# Create new arrays that share the same underlying data, but319# have 0 and 1 reversed for the first coordinate.320self.y_vertices_swapped = self.y_vertices[::-1, ...]321self.z_vertices_swapped = self.z_vertices[::-1, ...]322323self.slices = []324325def process_slice(self, unsigned int x, np.ndarray slice):326"""327Process a single slice of function evaluations at the specified x coordinate.328329EXAMPLES::330331sage: from sage.plot.plot3d.implicit_surface import MarchingCubesTriangles332sage: import numpy as np333sage: cube_marcher = MarchingCubesTriangles((-2, 2), (-2, 2), (-2, 2), 4, (10,)*3, smooth=False)334sage: f = lambda x, y, z: x^2 + y^2 + z^2335sage: slices = np.zeros((10, 10, 10), dtype=np.double)336sage: for x in reversed(xrange(0, 10)):337... for y in xrange(0, 10):338... for z in xrange(0, 10):339... slices[x, y, z] = f(*[a * (4 / 9) -2 for a in (x, y, z)])340... cube_marcher.process_slice(x, slices[x, :, :])341sage: faces = cube_marcher.finish()342sage: faces[0][0]343{'y': -1.1..., 'x': 1.5..., 'z': -0.5...}344345We render the isosurface using IndexFaceSet::346347sage: from sage.plot.plot3d.index_face_set import IndexFaceSet348sage: IndexFaceSet([tuple((p['x'], p['y'], p['z']) for p in face) for face in faces])349"""350# We go to a lot of effort to avoid repeating computations.351# (I hope that the related bookkeeping is actually faster352# than repeating the computations, but I haven't tested.)353354# We get evaluation points, one slice at a time. But we355# don't want to process slices of evaluation points;356# we want to process slices of cubes. (There is one fewer slice357# of cubes than there is of evaluation points.)358359# Without any caching, we would repeat this for each cube:360# Find which vertices are inside/outside the surface.361# For edges which cross the surface, interpolate to362# find the exact crossing location; these will be the vertices363# of triangles.364# Compute the color and surface normal for each of these vertices.365# Assemble the vertices into triangles, using the "marching cubes"366# tables.367368# The problem with this is that each vertex is actually shared369# among four neighbor cubes (except on the edge of the array),370# so we would be repeating work (potentially a lot of work, if371# the user-specified gradient or color routines are expensive)372# four times.373374# So we cache the information associated with each vertex.375# Let's call an edge an X, Y, or Z edge depending on which376# axis it is parallel to; and a vertex is an X, Y, or Z377# vertex depending on which kind of edge it lies on.378379# X vertices are shared between cubes that are all within a380# single cube slice. However, Y and Z vertices are shared381# between adjacent slices, so we need to keep those caches382# around. But we still need only two slices of vertex cache383# at a time (the two that are adjacent to the current cube slice).384385# To reduce memory allocations, we allocate these caches at the386# beginning of the run (they are ndarrays of objects). We387# set a VertexInfo object in the cache when we first compute it;388# then we look in the cache when we need it for the other three389# cubes.390391# The X vertex cache is a 2-D ndarray. The Y and Z vertex caches392# are 3-D, where the first (x) dimension is indexed by 0 or 1.393394# The Cython buffer interface (that we use for fast access to395# ndarray's) has to reinitialize itself separately in each396# function. For this reason, we use fewer, bigger functions;397# in particular, we don't call any functions per-cube that398# require buffer access.399400# To compute interpolated gradients using central differencing,401# we need up to 4 slices. Consider a gradient computed at402# an X vertex, between slices 1 and 2. This is interpolated403# between the gradient at slice 1 and the gradient at slice 2.404# Computing the gradient at slice 1 requires slices 0 and 2;405# computing the gradient at slice 2 requires slices 1 and 3.406407# This means we need to queue up slices and process them408# in a somewhat strange order. This function only does409# the queuing, and then passes all the real work off to410# other functions.411412self.slices = ([slice] + self.slices)[:4]413414if len(self.slices) >= 2:415self._update_yz_vertices(x+1,416self.slices[0],417self.slices[1],418self.slices[2] if len(self.slices) >= 3 else None)419420if len(self.slices) >= 3:421self._update_x_vertices(x+1,422self.slices[0],423self.slices[1],424self.slices[2],425self.slices[3] if len(self.slices) >= 4 else None)426self.process_cubes(self.slices[1], self.slices[2])427428if x == 0:429self._update_yz_vertices(x,430None,431self.slices[0],432self.slices[1])433434self._update_x_vertices(x,435None,436self.slices[0],437self.slices[1],438self.slices[2] if len(self.slices) >= 3 else None)439440self.process_cubes(self.slices[0], self.slices[1])441442cpdef _update_yz_vertices(self, int x, np.ndarray _prev, np.ndarray _cur, np.ndarray _next):443"""444TESTS::445446sage: from sage.plot.plot3d.implicit_surface import MarchingCubesTriangles447sage: import numpy as np448sage: cube_marcher = MarchingCubesTriangles((0, 1), (0, 1), (0, 1), 0, (2,)*3, smooth=False)449sage: prev_slice = next_slice = np.ones((2, 2), dtype=np.double)450sage: cur_slice = prev_slice.copy()451sage: cur_slice[0, 0] = -1 # Seed the slice data with an "interesting" value.452sage: cube_marcher._update_yz_vertices(1, prev_slice, cur_slice, next_slice)453sage: cube_marcher.z_vertices.tolist()454[[[<1.0, 0.0, 0.5>], [None]], [[None], [None]]]455sage: cube_marcher.y_vertices.tolist()456[[[<1.0, 0.5, 0.0>, None]], [[None, None]]]457sage: cube_marcher.x_vertices.any() # This shouldn't affect the X vertices.458False459"""460(self.y_vertices, self.y_vertices_swapped) = \461(self.y_vertices_swapped, self.y_vertices)462(self.z_vertices, self.z_vertices_swapped) = \463(self.z_vertices_swapped, self.z_vertices)464465cdef bint has_prev = (_prev is not None)466cdef bint has_next = (_next is not None)467468cdef np.ndarray[double, ndim=2] prev = _prev469cdef np.ndarray[double, ndim=2] cur = _cur470cdef np.ndarray[double, ndim=2] next = _next471472cdef np.ndarray[object, ndim=2] y_vertices = self.y_vertices[0,...]473cdef np.ndarray[object, ndim=2] z_vertices = self.z_vertices[0,...]474475cdef int ny = self.ny476cdef int nz = self.nz477478cdef int y479cdef int z480cdef VertexInfo v481cdef double frac482cdef bint has_nan483cdef point_c gradients[2]484cdef int i485for y from 0 <= y < ny - 1:486for z from 0 <= z < nz:487if marching_has_edge(cur[y,z], cur[y+1,z], self.contour, &frac, &has_nan):488v = mk_VertexInfo(x, y+frac, z, &self.eval_min, &self.eval_scale)489if self.region is not None:490if not self.in_region(v):491y_vertices[y,z] = None492continue493if self.smooth:494# We must compute a gradient.495if self.gradient is not None:496self.apply_point_func(&v.gradient, self.gradient, v)497else:498# Use central differencing.499for i from 0 <= i < 2:500self.get_gradient(&gradients[i],501x, y+i, z,502cur[y+i,z],503prev[y+i,z] if has_prev else 0,504next[y+i,z] if has_next else 0,505cur[y+i-1,z] if y+i>0 else 0,506cur[y+i+1,z] if y+i<ny-1 else 0,507cur[y+i,z-1] if z>0 else 0,508cur[y+i,z+1] if z<nz-1 else 0)509interpolate_point_c(&v.gradient, frac, gradients)510if self.vertex_color:511self.apply_point_func(&v.color, self.vertex_color, v)512y_vertices[y,z] = v513else:514y_vertices[y,z] = None515516# OK, that updated the Y vertices. Now we do almost exactly517# the same thing to update Z vertices.518for y from 0 <= y < ny:519for z from 0 <= z < nz - 1:520if marching_has_edge(cur[y,z], cur[y,z+1], self.contour, &frac, &has_nan):521v = mk_VertexInfo(x, y, z+frac, &self.eval_min, &self.eval_scale)522if self.region is not None:523if not self.in_region(v):524z_vertices[y,z] = None525continue526if self.smooth:527# We must compute a gradient.528if self.gradient is not None:529self.apply_point_func(&v.gradient, self.gradient, v)530else:531# Use central differencing.532for i from 0 <= i < 2:533self.get_gradient(&gradients[i],534x, y, z+i,535cur[y,z+i],536prev[y,z+i] if has_prev else 0,537next[y,z+i] if has_next else 0,538cur[y-1,z+i] if y>0 else 0,539cur[y+1,z+i] if y<ny-1 else 0,540cur[y,z+i-1] if z+i>0 else 0,541cur[y,z+i+1] if z+i<nz-1 else 0)542interpolate_point_c(&v.gradient, frac, gradients)543if self.vertex_color:544self.apply_point_func(&v.color, self.vertex_color, v)545z_vertices[y,z] = v546else:547z_vertices[y,z] = None548549cpdef _update_x_vertices(self, int x, np.ndarray _prev, np.ndarray _left, np.ndarray _right, np.ndarray _next):550"""551TESTS::552553sage: from sage.plot.plot3d.implicit_surface import MarchingCubesTriangles554sage: import numpy as np555sage: cube_marcher = MarchingCubesTriangles((0, 1), (0, 1), (0, 1), 0, (4, 2, 2), smooth=False)556sage: prev_slice = right_slice = next_slice = np.ones((2, 2), dtype=np.double)557sage: left_slice = prev_slice.copy()558sage: left_slice[1, 1] = -1559sage: cube_marcher._update_x_vertices(1, prev_slice, left_slice, right_slice, next_slice)560sage: cube_marcher.x_vertices.tolist()561[[None, None], [None, <1.5, 1.0, 1.0>]]562sage: cube_marcher.y_vertices.any() or cube_marcher.z_vertices.any() # This shouldn't affect the Y or Z vertices.563False564"""565cdef bint has_prev = (_prev is not None)566cdef bint has_next = (_next is not None)567568cdef np.ndarray[double, ndim=2] prev = _prev569cdef np.ndarray[double, ndim=2] left = _left570cdef np.ndarray[double, ndim=2] right = _right571cdef np.ndarray[double, ndim=2] next = _next572573cdef np.ndarray[object, ndim=2] x_vertices = self.x_vertices574575cdef int ny = self.ny576cdef int nz = self.nz577578cdef int y579cdef int z580581cdef VertexInfo v582cdef double frac583cdef bint has_nan584cdef point_c gradients[2]585for y from 0 <= y < ny:586for z from 0 <= z < nz:587if marching_has_edge(left[y,z], right[y,z], self.contour, &frac, &has_nan):588v = mk_VertexInfo(x+frac, y, z, &self.eval_min, &self.eval_scale)589if self.region is not None:590if not self.in_region(v):591x_vertices[y,z] = None592continue593if self.smooth:594# We must compute a gradient.595if self.gradient is not None:596self.apply_point_func(&v.gradient, self.gradient, v)597else:598# Use central differencing.599self.get_gradient(&gradients[0],600x, y, z,601left[y,z],602prev[y,z] if has_prev else 0,603right[y,z],604left[y-1,z] if y>0 else 0,605left[y+1,z] if y<ny-1 else 0,606left[y,z-1] if z>0 else 0,607left[y,z+1] if z<nz-1 else 0)608self.get_gradient(&gradients[1],609x, y, z,610right[y,z],611left[y,z],612next[y,z] if has_next else 0,613right[y-1,z] if y>0 else 0,614right[y+1,z] if y<ny-1 else 0,615right[y,z-1] if z>0 else 0,616right[y,z+1] if z<nz-1 else 0)617interpolate_point_c(&v.gradient, frac, gradients)618if self.vertex_color:619self.apply_point_func(&v.color, self.vertex_color, v)620x_vertices[y,z] = v621else:622x_vertices[y,z] = None623624cdef bint in_region(self, VertexInfo v):625return (self.region(v.eval_pt.x, v.eval_pt.y, v.eval_pt.z) > 0)626627cdef apply_point_func(self, point_c *pt, fn, VertexInfo v):628if isinstance(fn, tuple):629pt[0].x = fn[0](v.eval_pt.x, v.eval_pt.y, v.eval_pt.z)630pt[0].y = fn[1](v.eval_pt.x, v.eval_pt.y, v.eval_pt.z)631pt[0].z = fn[2](v.eval_pt.x, v.eval_pt.y, v.eval_pt.z)632else:633t = fn(v.eval_pt.x, v.eval_pt.y, v.eval_pt.z)634pt[0].x = t[0]635pt[0].y = t[1]636pt[0].z = t[2]637638cdef get_gradient(self,639point_c *g,640int x, int y, int z,641double center,642double lx, double ux,643double ly, double uy,644double lz, double uz):645# What a mess! It would be much nicer-looking code to pass slices646# in here and do the subscripting in here. Unfortunately,647# that would also be much slower, because we'd have to re-initialize648# the Cython buffer interface on each call.649650cdef double dx = ux - lx651cdef double dy = uy - ly652cdef double dz = uz - lz653654cdef double gx = dx * self.eval_scale_inv.x655cdef double gy = dy * self.eval_scale_inv.y656cdef double gz = dz * self.eval_scale_inv.z657658if x > 0 and x < self.nx - 1: gx *= 0.5659if y > 0 and y < self.ny - 1: gy *= 0.5660if z > 0 and z < self.nz - 1: gz *= 0.5661662g[0].x = gx663g[0].y = gy664g[0].z = gz665666cpdef process_cubes(self, np.ndarray _left, np.ndarray _right):667"""668TESTS::669670sage: from sage.plot.plot3d.implicit_surface import MarchingCubesTriangles671sage: import numpy as np672sage: cube_marcher = MarchingCubesTriangles((0, 1), (0, 1), (0, 1), 0, (3, 2, 2), smooth=False)673sage: slices = [np.ones((2, 2), dtype=np.double) for i in xrange(0, 3)]674sage: slices[0][1, 1] = -1675sage: cube_marcher._update_yz_vertices(0, None, slices[0], slices[1])676sage: cube_marcher._update_x_vertices(0, None, slices[0], slices[1], slices[2])677sage: cube_marcher.process_cubes(slices[0], slices[1])678sage: cube_marcher.finish()679[({'y': 1.0, 'x': 0.0, 'z': 0.5}, {'y': 1.0, 'x': 0.25, 'z': 1.0}, {'y': 0.5, 'x': 0.0, 'z': 1.0})]680"""681cdef np.ndarray[double, ndim=2] left = _left682cdef np.ndarray[double, ndim=2] right = _right683684cdef np.ndarray[object, ndim=2] x_vertices = self.x_vertices685cdef np.ndarray[object, ndim=3] y_vertices = self.y_vertices686cdef np.ndarray[object, ndim=3] z_vertices = self.z_vertices687688cdef int ny = self.ny689cdef int nz = self.nz690691cdef int y692cdef int z693694# based on generateSurfaceData in MarchingCubes.java695cdef int insideMask696697# Cool ASCII art from MarchingCubes.java:698# * Y699# * 4 --------4--------- 5700# * /| /|701# * / | / |702# * / | / |703# * 7 8 5 |704# * / | / 9705# * / | / |706# * 7 --------6--------- 6 |707# * | | | |708# * | 0 ---------0--|----- 1 X709# * | / | /710# * 11 / 10 /711# * | 3 | 1712# * | / | /713# * | / | /714# * 3 ---------2-------- 2715# * Z716717# We see from the above that vertices are labeled 0 to 7, and718# edges are labeled 0 to 11.719720cdef list all_vertex_info = [None] * 12721cdef tuple my_triangles722723cdef int i724725for y from 0 <= y < ny-1:726for z from 0 <= z < nz-1:727# For each vertex (0 to 7), set the corresponding bit728# of insideMask iff the vertex is inside the surface.729insideMask = 0730insideMask |= marching_is_inside(left[y,z], self.contour)<<0731insideMask |= marching_is_inside(right[y,z], self.contour)<<1732insideMask |= marching_is_inside(right[y,z+1], self.contour)<<2733insideMask |= marching_is_inside(left[y,z+1], self.contour)<<3734insideMask |= marching_is_inside(left[y+1,z], self.contour)<<4735insideMask |= marching_is_inside(right[y+1,z], self.contour)<<5736insideMask |= marching_is_inside(right[y+1,z+1], self.contour)<<6737insideMask |= marching_is_inside(left[y+1,z+1], self.contour)<<7738739if insideMask == 0 or insideMask == 255: continue740741# OK, we have a cube on the surface. Copy all of the vertex742# info into an array for easier reference.743744all_vertex_info[0] = x_vertices[y,z]745all_vertex_info[1] = z_vertices[1,y,z]746all_vertex_info[2] = x_vertices[y,z+1]747all_vertex_info[3] = z_vertices[0,y,z]748all_vertex_info[4] = x_vertices[y+1,z]749all_vertex_info[5] = z_vertices[1,y+1,z]750all_vertex_info[6] = x_vertices[y+1,z+1]751all_vertex_info[7] = z_vertices[0,y+1,z]752all_vertex_info[8] = y_vertices[0,y,z]753all_vertex_info[9] = y_vertices[1,y,z]754all_vertex_info[10]= y_vertices[1,y,z+1]755all_vertex_info[11]= y_vertices[0,y,z+1]756757my_triangles = triangle_table2[insideMask]758759for i in range(0, len(my_triangles), 4):760# In wireframe mode, my_triangles[i+3] specifies761# whether or not to draw the corresponding edge.762# See MarchingCubes.java for details.763self.add_triangle(all_vertex_info[my_triangles[i]],764all_vertex_info[my_triangles[i+1]],765all_vertex_info[my_triangles[i+2]])766767cpdef add_triangle(self, VertexInfo v1, VertexInfo v2, VertexInfo v3):768"""769Called when a new triangle is generated by the marching cubes algorithm770to update the results array.771772TESTS::773774sage: from sage.plot.plot3d.implicit_surface import MarchingCubesTriangles, VertexInfo775sage: cube_marcher = MarchingCubesTriangles((0, 1), (0, 1), (0, 1), 0, (10,)*3, smooth=False)776sage: cube_marcher.add_triangle(VertexInfo(), VertexInfo(), VertexInfo())777sage: cube_marcher.finish()778[({'y': 0.0, 'x': 0.0, 'z': 0.0}, {'y': 0.0, 'x': 0.0, 'z': 0.0}, {'y': 0.0, 'x': 0.0, 'z': 0.0})]779"""780781if v1 is None or v2 is None or v3 is None:782# This happens if there is a NaN nearby, or if a hole was specified here.783return784785cdef:786point_c v1_ev_pt, v2_ev_pt, v3_ev_pt787point_c n1_ev_vec, n2_ev_vec, n3_ev_vec788789if self.transform is not None:790self.transform.transform_point_c(&v1_ev_pt, v1.eval_pt)791self.transform.transform_point_c(&v2_ev_pt, v2.eval_pt)792self.transform.transform_point_c(&v3_ev_pt, v3.eval_pt)793else:794v1_ev_pt = v1.eval_pt795v2_ev_pt = v2.eval_pt796v3_ev_pt = v3.eval_pt797face = (v1_ev_pt, v2_ev_pt, v3_ev_pt)798799if self.smooth:800# XXX I believe this is wrong for non-uniform transforms801if self.transform is not None:802self.transform.transform_vector_c(&n1_ev_vec, v1.gradient)803self.transform.transform_vector_c(&n2_ev_vec, v2.gradient)804self.transform.transform_vector_c(&n3_ev_vec, v3.gradient)805else:806n1_ev_vec = v1.gradient807n2_ev_vec = v2.gradient808n3_ev_vec = v3.gradient809face += (n1_ev_vec, n2_ev_vec, n3_ev_vec)810811self.results.append(face)812813cpdef render_implicit(f, xrange, yrange, zrange, plot_points, cube_marchers):814"""815INPUT:816817- ``f`` - a (fast!) callable function818819- ``xrange`` - a 2-tuple (x_min, x_max)820821- ``yrange`` - a 2-tuple (y_min, y_may)822823- ``zrange`` - a 2-tuple (z_min, z_maz)824825- ``plot_points`` - a triple of integers indicating the number of826function evaluations in each direction.827828- ``cube_marchers`` - a list of cube marchers, one for each contour.829830831OUTPUT:832833A representation of the isosurface, in the format specified by the individual834cube marchers.835836TESTS::837838sage: from sage.plot.plot3d.implicit_surface import render_implicit, MarchingCubesTriangles839sage: plot_points, f = (40,)*3, lambda x, y, z: x + y + z840sage: cube_marcher = MarchingCubesTriangles((0, 1), (0, 1), (0, 1), 1, (10,)*3)841sage: results = render_implicit(lambda x, y, z: x + y + z, \842... (0, 1), (0, 1), (0, 1), (10,)*3, [cube_marcher])843sage: results[0][0]844{'y': 0.0, 'x': 1.0, 'z': 0.0}845"""846847cdef int nx = plot_points[0]848cdef int ny = plot_points[1]849cdef int nz = plot_points[2]850851cdef point_c eval_min, eval_scale852853eval_min.x = xrange[0]854eval_scale.x = (xrange[1] - xrange[0]) / (nx - 1)855eval_min.y = yrange[0]856eval_scale.y = (yrange[1] - yrange[0]) / (ny - 1)857eval_min.z = zrange[0]858eval_scale.z = (zrange[1] - zrange[0]) / (nz - 1)859860# A possible future extension would be to allow passing in "f"861# as a numpy ndarray. If that were done, we could slightly modify862# the following code to just pass slices of f to the renderers863# (no copying of the data would be required).864865# The current marching cube renderers need only at most four slices at866# a time.867868cdef np.ndarray[double, ndim=3] data = np.zeros((4, ny, nz), dtype=np.double)869cdef np.ndarray[double, ndim=2] slice870871cdef unsigned int x, y, z872cdef int n873874cdef double eval_x, eval_y, eval_z875876cdef MarchingCubes marcher877878for n from 0 <= n < nx:879x = nx-1-n880eval_x = eval_min.x + eval_scale.x * x881slice = data[n % 4, :, :]882for y from 0 <= y < ny:883eval_y = eval_min.y + eval_scale.y * y884for z from 0 <= z < nz:885eval_z = eval_min.z + eval_scale.z * z886slice[y, z] = f(eval_x, eval_y, eval_z)887888for marcher in cube_marchers:889marcher.process_slice(x, slice)890891results = []892893for marcher in cube_marchers:894results.extend(marcher.finish())895896return results897898cdef class ImplicitSurface(IndexFaceSet):899cdef readonly object f900cdef readonly object vars901cdef readonly tuple xrange902cdef readonly tuple yrange903cdef readonly tuple zrange904cdef readonly list contours905cdef readonly object region906cdef readonly bint smooth907cdef readonly object gradient908cdef readonly object vertex_color909cdef readonly tuple plot_points910911def __init__(self, f, xrange, yrange, zrange,912contour=0, plot_points="automatic",913region=None, smooth=True, gradient=None, vertex_color=None,914**kwds):915"""916TESTS::917918sage: from sage.plot.plot3d.implicit_surface import ImplicitSurface919sage: var('x,y,z')920(x, y, z)921sage: G = ImplicitSurface(x^2 + y^2 + z^2, (x,-2, 2), (y,-2, 2), (z,-2, 2), contour=4)922sage: show(G)923"""924IndexFaceSet.__init__(self, [], [], **kwds)925from sage.ext.fast_eval import fast_float926927orig_f = f928self.f, ranges, self.vars = setup_for_eval_on_grid(f, [xrange, yrange, zrange], return_vars=True)929self.xrange = ranges[0][:2]930self.yrange = ranges[1][:2]931self.zrange = ranges[2][:2]932if isinstance(contour, (list, tuple)):933contours = contour934else:935contours = [contour]936self.contours = [float(c) for c in contours]937if region is not None:938self.region = fast_float(region, *self.vars)939940# Comments from Carl Witty, who first wrote this some of this code941# See Trac 9483942# When I first wrote the code, I had the idea to create a943# direct-to-tachyon backend that would use vertex normals944# to create much nicer-looking plots with smaller numbers945# of plot_points, and laid the groundwork for this backend946# with the gradient and smooth arguments. But I abandoned the947# code without writing this backend (and leaving many other parts948# of the code unfinished).949# When William Cauchois took over and finished the code (thank you,950# William!), he only wrote an IndexFaceSet backend, that can't951# (currently) make use of vertex normals. So the gradient code is952# currently useless.953# But it's still open for somebody to write a direct-to-tachyon backend,954# or to extend IndexFaceSet to support vertex normals.955956# Since IndexFaceSet doesn't even support smooth shading, we overwrite957# the passed-in smooth parameter.958smooth=False959960961self.smooth = smooth962if smooth and gradient is None:963try:964gradient = (orig_f.diff(self.vars[0]),965orig_f.diff(self.vars[1]),966orig_f.diff(self.vars[2]))967except:968# Would be nice to have more nuanced error handling here.969970# If anything goes wrong, we'll just use central differencing.971pass972if gradient is not None:973self.gradient = fast_float(gradient, *self.vars)974if vertex_color is not None:975self.vertex_color = fast_float(vertex_color, *self.vars)976if plot_points == "automatic":977plot_points = DEFAULT_PLOT_POINTS978my_plot_points = []979for i in range(3):980if isinstance(plot_points, (list, tuple)):981n = int(plot_points[i])982else:983n = int(plot_points)984if n < 2:985raise ValueError986my_plot_points.append(n)987self.plot_points = tuple(my_plot_points)988989def bounding_box(self):990"""991Return a bounding box for the ``ImplicitSurface``, as a tuple of two9923-dimensional points.993994EXAMPLES:995996Note that the bounding box corresponds exactly to the x-, y-, and z- range::997998sage: from sage.plot.plot3d.implicit_surface import ImplicitSurface999sage: G = ImplicitSurface(0, (0, 1), (0, 1), (0, 1))1000sage: G.bounding_box()1001((0.0, 0.0, 0.0), (1.0, 1.0, 1.0))1002"""1003return ((self.xrange[0], self.yrange[0], self.zrange[0]),1004(self.xrange[1], self.yrange[1], self.zrange[1]))10051006def obj_repr(self, render_params):1007"""1008Return a representation of this object in the .obj format.10091010TESTS::10111012We graph a simple plane::10131014sage: from sage.plot.plot3d.implicit_surface import ImplicitSurface1015sage: var('x,y,z')1016(x, y, z)1017sage: G = ImplicitSurface(x + y + z, (x,-1, 1), (y,-1, 1), (z,-1, 1))1018sage: obj = G.obj_repr(G.default_render_params())1019sage: vertices = obj[2]10201021The number of vertices in the OBJ representation should equal the number1022of vertices in the face set::10231024sage: len(vertices) == len(G.vertex_list())1025True10261027The vertices in the OBJ representation should also be approximately equal1028to the vertices in the face set -- the small error is due to rounding1029which occurs during output (we test only the first 20 points for the1030sake of speed)::10311032sage: def points_equal(a, b, epsilon=(1e-5)):1033... return all(abs(x0-x1) < epsilon for x0, x1 in zip(a, b))1034sage: list = []1035sage: assert len(vertices) >= 20 # I should hope so, we're rendering at the default resolution!1036sage: for vertex, surf_vertex in zip(vertices, G.vertex_list())[0:20]:1037... list.append(points_equal(map(float, vertex.split(' ')[1:]), surf_vertex))1038sage: all(list)1039True1040"""1041self.triangulate()1042return IndexFaceSet.obj_repr(self, render_params)10431044def tachyon_repr(self, render_params):1045"""1046Return a representation of this object suitable for use with the Tachyon1047renderer.10481049TESTS::10501051sage: from sage.plot.plot3d.implicit_surface import ImplicitSurface1052sage: var('x,y,z')1053(x, y, z)1054sage: G = ImplicitSurface(x + y + z, (x,-1, 1), (y,-1, 1), (z,-1, 1))1055sage: G.tachyon_repr(G.default_render_params())[0].startswith('TRI')1056True1057"""1058self.triangulate()1059return IndexFaceSet.tachyon_repr(self, render_params)10601061def jmol_repr(self, render_params):1062"""1063Return a representation of this object suitable for use with the Jmol1064renderer.10651066TESTS::10671068sage: from sage.plot.plot3d.implicit_surface import ImplicitSurface1069sage: var('x,y,z')1070(x, y, z)1071sage: G = ImplicitSurface(x + y + z, (x,-1, 1), (y,-1, 1), (z,-1, 1))1072sage: show(G, viewer='jmol')1073"""1074self.triangulate()1075return IndexFaceSet.jmol_repr(self, render_params)10761077def json_repr(self, render_params):1078"""1079Return a representation of this object in JavaScript Object Notation (JSON).10801081TESTS::10821083sage: from sage.plot.plot3d.implicit_surface import ImplicitSurface1084sage: var('x,y,z')1085(x, y, z)1086sage: G = ImplicitSurface(x + y + z, (x,-1, 1), (y,-1, 1), (z,-1, 1))1087sage: G.json_repr(G.default_render_params())[0].startswith('{vertices:')1088True1089"""1090self.triangulate()1091return IndexFaceSet.json_repr(self, render_params)10921093def triangulate(self, force=False):1094"""1095The IndexFaceSet will be empty until you call this method, which generates1096the faces and vertices according to the parameters specified in the1097constructor for ImplicitSurface. Note that if you call this method more1098than once, subsequent invocations will have no effect (this is an1099optimization to avoid repeated work) unless you specify force=True in the1100keywords.11011102EXAMPLES::11031104sage: from sage.plot.plot3d.implicit_surface import ImplicitSurface1105sage: var('x,y,z')1106(x, y, z)1107sage: G = ImplicitSurface(x + y + z, (x,-1, 1), (y,-1, 1), (z,-1, 1))1108sage: len(G.vertex_list()), len(G.face_list())1109(0, 0)1110sage: G.triangulate()1111sage: len(G.vertex_list()) > 0, len(G.face_list()) > 01112(True, True)1113sage: G.show() # This should be fast, since the mesh is already triangulated.1114"""1115if self.fcount != 0 and not force:1116# The mesh is already triangulated1117return11181119options = dict(xrange=self.xrange, yrange=self.yrange, zrange=self.zrange,1120region=self.region, smooth=self.smooth,1121gradient=self.gradient,1122vertex_color=self.vertex_color,1123plot_points=self.plot_points)1124cube_marchers = [MarchingCubesTriangles(contour=x, **options) for x in self.contours]1125results = render_implicit(self.f, self.xrange, self.yrange, self.zrange,1126self.plot_points, cube_marchers)1127cdef:1128face_c* dest_face1129point_c* dest_vertex1130int fcount = len(results)11311132self.realloc(fcount * 3, fcount, fcount * 3)1133for i from 0 <= i < fcount:1134dest_face = &self._faces[i]1135src_face = results[i]11361137dest_face.n = 31138dest_face.vertices = &self.face_indices[3 * i]11391140for j from 0 <= j < 3:1141dest_face.vertices[j] = (3 * i) + j1142dest_vertex = &self.vs[(3 * i) + j]1143dest_vertex.x = src_face[j]['x']1144dest_vertex.y = src_face[j]['y']1145dest_vertex.z = src_face[j]['z']11461147self.vcount = fcount * 31148self.fcount = fcount1149self.icount = fcount * 311501151# Data table (courtesy of MarchingCubes.java)1152triangle_table2 = ( None,1153( 0, 8, 3, 7 ),1154( 0, 1, 9, 7 ), ( 1, 8, 3, 6, 9, 8, 1, 5 ), ( 1, 2, 10, 7 ),1155( 0, 8, 3, 7, 1, 2, 10, 7 ), ( 9, 2, 10, 6, 0, 2, 9, 5 ),1156( 2, 8, 3, 6, 2, 10, 8, 1, 10, 9, 8, 3 ), ( 3, 11, 2, 7 ),1157( 0, 11, 2, 6, 8, 11, 0, 5 ), ( 1, 9, 0, 7, 2, 3, 11, 7 ),1158( 1, 11, 2, 6, 1, 9, 11, 1, 9, 8, 11, 3 ), ( 3, 10, 1, 6, 11, 10, 3, 5 ),1159( 0, 10, 1, 6, 0, 8, 10, 1, 8, 11, 10, 3 ),1160( 3, 9, 0, 6, 3, 11, 9, 1, 11, 10, 9, 3 ), ( 9, 8, 10, 5, 10, 8, 11, 6 ),1161( 4, 7, 8, 7 ), ( 4, 3, 0, 6, 7, 3, 4, 5 ), ( 0, 1, 9, 7, 8, 4, 7, 7 ),1162( 4, 1, 9, 6, 4, 7, 1, 1, 7, 3, 1, 3 ), ( 1, 2, 10, 7, 8, 4, 7, 7 ),1163( 3, 4, 7, 6, 3, 0, 4, 3, 1, 2, 10, 7 ),1164( 9, 2, 10, 6, 9, 0, 2, 3, 8, 4, 7, 7 ),1165( 2, 10, 9, 3, 2, 9, 7, 0, 2, 7, 3, 6, 7, 9, 4, 6 ),1166( 8, 4, 7, 7, 3, 11, 2, 7 ), ( 11, 4, 7, 6, 11, 2, 4, 1, 2, 0, 4, 3 ),1167( 9, 0, 1, 7, 8, 4, 7, 7, 2, 3, 11, 7 ),1168( 4, 7, 11, 3, 9, 4, 11, 1, 9, 11, 2, 2, 9, 2, 1, 6 ),1169( 3, 10, 1, 6, 3, 11, 10, 3, 7, 8, 4, 7 ),1170( 1, 11, 10, 6, 1, 4, 11, 0, 1, 0, 4, 3, 7, 11, 4, 5 ),1171( 4, 7, 8, 7, 9, 0, 11, 1, 9, 11, 10, 6, 11, 0, 3, 6 ),1172( 4, 7, 11, 3, 4, 11, 9, 4, 9, 11, 10, 6 ), ( 9, 5, 4, 7 ),1173( 9, 5, 4, 7, 0, 8, 3, 7 ), ( 0, 5, 4, 6, 1, 5, 0, 5 ),1174( 8, 5, 4, 6, 8, 3, 5, 1, 3, 1, 5, 3 ), ( 1, 2, 10, 7, 9, 5, 4, 7 ),1175( 3, 0, 8, 7, 1, 2, 10, 7, 4, 9, 5, 7 ),1176( 5, 2, 10, 6, 5, 4, 2, 1, 4, 0, 2, 3 ),1177( 2, 10, 5, 3, 3, 2, 5, 1, 3, 5, 4, 2, 3, 4, 8, 6 ),1178( 9, 5, 4, 7, 2, 3, 11, 7 ), ( 0, 11, 2, 6, 0, 8, 11, 3, 4, 9, 5, 7 ),1179( 0, 5, 4, 6, 0, 1, 5, 3, 2, 3, 11, 7 ),1180( 2, 1, 5, 3, 2, 5, 8, 0, 2, 8, 11, 6, 4, 8, 5, 5 ),1181( 10, 3, 11, 6, 10, 1, 3, 3, 9, 5, 4, 7 ),1182( 4, 9, 5, 7, 0, 8, 1, 5, 8, 10, 1, 2, 8, 11, 10, 3 ),1183( 5, 4, 0, 3, 5, 0, 11, 0, 5, 11, 10, 6, 11, 0, 3, 6 ),1184( 5, 4, 8, 3, 5, 8, 10, 4, 10, 8, 11, 6 ), ( 9, 7, 8, 6, 5, 7, 9, 5 ),1185( 9, 3, 0, 6, 9, 5, 3, 1, 5, 7, 3, 3 ),1186( 0, 7, 8, 6, 0, 1, 7, 1, 1, 5, 7, 3 ), ( 1, 5, 3, 5, 3, 5, 7, 6 ),1187( 9, 7, 8, 6, 9, 5, 7, 3, 10, 1, 2, 7 ),1188( 10, 1, 2, 7, 9, 5, 0, 5, 5, 3, 0, 2, 5, 7, 3, 3 ),1189( 8, 0, 2, 3, 8, 2, 5, 0, 8, 5, 7, 6, 10, 5, 2, 5 ),1190( 2, 10, 5, 3, 2, 5, 3, 4, 3, 5, 7, 6 ),1191( 7, 9, 5, 6, 7, 8, 9, 3, 3, 11, 2, 7 ),1192( 9, 5, 7, 3, 9, 7, 2, 0, 9, 2, 0, 6, 2, 7, 11, 6 ),1193( 2, 3, 11, 7, 0, 1, 8, 5, 1, 7, 8, 2, 1, 5, 7, 3 ),1194( 11, 2, 1, 3, 11, 1, 7, 4, 7, 1, 5, 6 ),1195( 9, 5, 8, 5, 8, 5, 7, 6, 10, 1, 3, 3, 10, 3, 11, 6 ),1196( 5, 7, 0, 1, 5, 0, 9, 6, 7, 11, 0, 1, 1, 0, 10, 5, 11, 10, 0, 1 ),1197( 11, 10, 0, 1, 11, 0, 3, 6, 10, 5, 0, 1, 8, 0, 7, 5, 5, 7, 0, 1 ),1198( 11, 10, 5, 3, 7, 11, 5, 5 ), ( 10, 6, 5, 7 ),1199( 0, 8, 3, 7, 5, 10, 6, 7 ), ( 9, 0, 1, 7, 5, 10, 6, 7 ),1200( 1, 8, 3, 6, 1, 9, 8, 3, 5, 10, 6, 7 ), ( 1, 6, 5, 6, 2, 6, 1, 5 ),1201( 1, 6, 5, 6, 1, 2, 6, 3, 3, 0, 8, 7 ),1202( 9, 6, 5, 6, 9, 0, 6, 1, 0, 2, 6, 3 ),1203( 5, 9, 8, 3, 5, 8, 2, 0, 5, 2, 6, 6, 3, 2, 8, 5 ),1204( 2, 3, 11, 7, 10, 6, 5, 7 ), ( 11, 0, 8, 6, 11, 2, 0, 3, 10, 6, 5, 7 ),1205( 0, 1, 9, 7, 2, 3, 11, 7, 5, 10, 6, 7 ),1206( 5, 10, 6, 7, 1, 9, 2, 5, 9, 11, 2, 2, 9, 8, 11, 3 ),1207( 6, 3, 11, 6, 6, 5, 3, 1, 5, 1, 3, 3 ),1208( 0, 8, 11, 3, 0, 11, 5, 0, 0, 5, 1, 6, 5, 11, 6, 6 ),1209( 3, 11, 6, 3, 0, 3, 6, 1, 0, 6, 5, 2, 0, 5, 9, 6 ),1210( 6, 5, 9, 3, 6, 9, 11, 4, 11, 9, 8, 6 ), ( 5, 10, 6, 7, 4, 7, 8, 7 ),1211( 4, 3, 0, 6, 4, 7, 3, 3, 6, 5, 10, 7 ),1212( 1, 9, 0, 7, 5, 10, 6, 7, 8, 4, 7, 7 ),1213( 10, 6, 5, 7, 1, 9, 7, 1, 1, 7, 3, 6, 7, 9, 4, 6 ),1214( 6, 1, 2, 6, 6, 5, 1, 3, 4, 7, 8, 7 ),1215( 1, 2, 5, 5, 5, 2, 6, 6, 3, 0, 4, 3, 3, 4, 7, 6 ),1216( 8, 4, 7, 7, 9, 0, 5, 5, 0, 6, 5, 2, 0, 2, 6, 3 ),1217( 7, 3, 9, 1, 7, 9, 4, 6, 3, 2, 9, 1, 5, 9, 6, 5, 2, 6, 9, 1 ),1218( 3, 11, 2, 7, 7, 8, 4, 7, 10, 6, 5, 7 ),1219( 5, 10, 6, 7, 4, 7, 2, 1, 4, 2, 0, 6, 2, 7, 11, 6 ),1220( 0, 1, 9, 7, 4, 7, 8, 7, 2, 3, 11, 7, 5, 10, 6, 7 ),1221( 9, 2, 1, 6, 9, 11, 2, 2, 9, 4, 11, 1, 7, 11, 4, 5, 5, 10, 6, 7 ),1222( 8, 4, 7, 7, 3, 11, 5, 1, 3, 5, 1, 6, 5, 11, 6, 6 ),1223( 5, 1, 11, 1, 5, 11, 6, 6, 1, 0, 11, 1, 7, 11, 4, 5, 0, 4, 11, 1 ),1224( 0, 5, 9, 6, 0, 6, 5, 2, 0, 3, 6, 1, 11, 6, 3, 5, 8, 4, 7, 7 ),1225( 6, 5, 9, 3, 6, 9, 11, 4, 4, 7, 9, 5, 7, 11, 9, 1 ),1226( 10, 4, 9, 6, 6, 4, 10, 5 ), ( 4, 10, 6, 6, 4, 9, 10, 3, 0, 8, 3, 7 ),1227( 10, 0, 1, 6, 10, 6, 0, 1, 6, 4, 0, 3 ),1228( 8, 3, 1, 3, 8, 1, 6, 0, 8, 6, 4, 6, 6, 1, 10, 6 ),1229( 1, 4, 9, 6, 1, 2, 4, 1, 2, 6, 4, 3 ),1230( 3, 0, 8, 7, 1, 2, 9, 5, 2, 4, 9, 2, 2, 6, 4, 3 ),1231( 0, 2, 4, 5, 4, 2, 6, 6 ), ( 8, 3, 2, 3, 8, 2, 4, 4, 4, 2, 6, 6 ),1232( 10, 4, 9, 6, 10, 6, 4, 3, 11, 2, 3, 7 ),1233( 0, 8, 2, 5, 2, 8, 11, 6, 4, 9, 10, 3, 4, 10, 6, 6 ),1234( 3, 11, 2, 7, 0, 1, 6, 1, 0, 6, 4, 6, 6, 1, 10, 6 ),1235( 6, 4, 1, 1, 6, 1, 10, 6, 4, 8, 1, 1, 2, 1, 11, 5, 8, 11, 1, 1 ),1236( 9, 6, 4, 6, 9, 3, 6, 0, 9, 1, 3, 3, 11, 6, 3, 5 ),1237( 8, 11, 1, 1, 8, 1, 0, 6, 11, 6, 1, 1, 9, 1, 4, 5, 6, 4, 1, 1 ),1238( 3, 11, 6, 3, 3, 6, 0, 4, 0, 6, 4, 6 ), ( 6, 4, 8, 3, 11, 6, 8, 5 ),1239( 7, 10, 6, 6, 7, 8, 10, 1, 8, 9, 10, 3 ),1240( 0, 7, 3, 6, 0, 10, 7, 0, 0, 9, 10, 3, 6, 7, 10, 5 ),1241( 10, 6, 7, 3, 1, 10, 7, 1, 1, 7, 8, 2, 1, 8, 0, 6 ),1242( 10, 6, 7, 3, 10, 7, 1, 4, 1, 7, 3, 6 ),1243( 1, 2, 6, 3, 1, 6, 8, 0, 1, 8, 9, 6, 8, 6, 7, 6 ),1244( 2, 6, 9, 1, 2, 9, 1, 6, 6, 7, 9, 1, 0, 9, 3, 5, 7, 3, 9, 1 ),1245( 7, 8, 0, 3, 7, 0, 6, 4, 6, 0, 2, 6 ), ( 7, 3, 2, 3, 6, 7, 2, 5 ),1246( 2, 3, 11, 7, 10, 6, 8, 1, 10, 8, 9, 6, 8, 6, 7, 6 ),1247( 2, 0, 7, 1, 2, 7, 11, 6, 0, 9, 7, 1, 6, 7, 10, 5, 9, 10, 7, 1 ),1248( 1, 8, 0, 6, 1, 7, 8, 2, 1, 10, 7, 1, 6, 7, 10, 5, 2, 3, 11, 7 ),1249( 11, 2, 1, 3, 11, 1, 7, 4, 10, 6, 1, 5, 6, 7, 1, 1 ),1250( 8, 9, 6, 1, 8, 6, 7, 6, 9, 1, 6, 1, 11, 6, 3, 5, 1, 3, 6, 1 ),1251( 0, 9, 1, 7, 11, 6, 7, 7 ),1252( 7, 8, 0, 3, 7, 0, 6, 4, 3, 11, 0, 5, 11, 6, 0, 1 ), ( 7, 11, 6, 7 ),1253( 7, 6, 11, 7 ), ( 3, 0, 8, 7, 11, 7, 6, 7 ),1254( 0, 1, 9, 7, 11, 7, 6, 7 ), ( 8, 1, 9, 6, 8, 3, 1, 3, 11, 7, 6, 7 ),1255( 10, 1, 2, 7, 6, 11, 7, 7 ), ( 1, 2, 10, 7, 3, 0, 8, 7, 6, 11, 7, 7 ),1256( 2, 9, 0, 6, 2, 10, 9, 3, 6, 11, 7, 7 ),1257( 6, 11, 7, 7, 2, 10, 3, 5, 10, 8, 3, 2, 10, 9, 8, 3 ),1258( 7, 2, 3, 6, 6, 2, 7, 5 ), ( 7, 0, 8, 6, 7, 6, 0, 1, 6, 2, 0, 3 ),1259( 2, 7, 6, 6, 2, 3, 7, 3, 0, 1, 9, 7 ),1260( 1, 6, 2, 6, 1, 8, 6, 0, 1, 9, 8, 3, 8, 7, 6, 3 ),1261( 10, 7, 6, 6, 10, 1, 7, 1, 1, 3, 7, 3 ),1262( 10, 7, 6, 6, 1, 7, 10, 4, 1, 8, 7, 2, 1, 0, 8, 3 ),1263( 0, 3, 7, 3, 0, 7, 10, 0, 0, 10, 9, 6, 6, 10, 7, 5 ),1264( 7, 6, 10, 3, 7, 10, 8, 4, 8, 10, 9, 6 ), ( 6, 8, 4, 6, 11, 8, 6, 5 ),1265( 3, 6, 11, 6, 3, 0, 6, 1, 0, 4, 6, 3 ),1266( 8, 6, 11, 6, 8, 4, 6, 3, 9, 0, 1, 7 ),1267( 9, 4, 6, 3, 9, 6, 3, 0, 9, 3, 1, 6, 11, 3, 6, 5 ),1268( 6, 8, 4, 6, 6, 11, 8, 3, 2, 10, 1, 7 ),1269( 1, 2, 10, 7, 3, 0, 11, 5, 0, 6, 11, 2, 0, 4, 6, 3 ),1270( 4, 11, 8, 6, 4, 6, 11, 3, 0, 2, 9, 5, 2, 10, 9, 3 ),1271( 10, 9, 3, 1, 10, 3, 2, 6, 9, 4, 3, 1, 11, 3, 6, 5, 4, 6, 3, 1 ),1272( 8, 2, 3, 6, 8, 4, 2, 1, 4, 6, 2, 3 ), ( 0, 4, 2, 5, 4, 6, 2, 3 ),1273( 1, 9, 0, 7, 2, 3, 4, 1, 2, 4, 6, 6, 4, 3, 8, 6 ),1274( 1, 9, 4, 3, 1, 4, 2, 4, 2, 4, 6, 6 ),1275( 8, 1, 3, 6, 8, 6, 1, 0, 8, 4, 6, 3, 6, 10, 1, 3 ),1276( 10, 1, 0, 3, 10, 0, 6, 4, 6, 0, 4, 6 ),1277( 4, 6, 3, 1, 4, 3, 8, 6, 6, 10, 3, 1, 0, 3, 9, 5, 10, 9, 3, 1 ),1278( 10, 9, 4, 3, 6, 10, 4, 5 ), ( 4, 9, 5, 7, 7, 6, 11, 7 ),1279( 0, 8, 3, 7, 4, 9, 5, 7, 11, 7, 6, 7 ),1280( 5, 0, 1, 6, 5, 4, 0, 3, 7, 6, 11, 7 ),1281( 11, 7, 6, 7, 8, 3, 4, 5, 3, 5, 4, 2, 3, 1, 5, 3 ),1282( 9, 5, 4, 7, 10, 1, 2, 7, 7, 6, 11, 7 ),1283( 6, 11, 7, 7, 1, 2, 10, 7, 0, 8, 3, 7, 4, 9, 5, 7 ),1284( 7, 6, 11, 7, 5, 4, 10, 5, 4, 2, 10, 2, 4, 0, 2, 3 ),1285( 3, 4, 8, 6, 3, 5, 4, 2, 3, 2, 5, 1, 10, 5, 2, 5, 11, 7, 6, 7 ),1286( 7, 2, 3, 6, 7, 6, 2, 3, 5, 4, 9, 7 ),1287( 9, 5, 4, 7, 0, 8, 6, 1, 0, 6, 2, 6, 6, 8, 7, 6 ),1288( 3, 6, 2, 6, 3, 7, 6, 3, 1, 5, 0, 5, 5, 4, 0, 3 ),1289( 6, 2, 8, 1, 6, 8, 7, 6, 2, 1, 8, 1, 4, 8, 5, 5, 1, 5, 8, 1 ),1290( 9, 5, 4, 7, 10, 1, 6, 5, 1, 7, 6, 2, 1, 3, 7, 3 ),1291( 1, 6, 10, 6, 1, 7, 6, 2, 1, 0, 7, 1, 8, 7, 0, 5, 9, 5, 4, 7 ),1292( 4, 0, 10, 1, 4, 10, 5, 6, 0, 3, 10, 1, 6, 10, 7, 5, 3, 7, 10, 1 ),1293( 7, 6, 10, 3, 7, 10, 8, 4, 5, 4, 10, 5, 4, 8, 10, 1 ),1294( 6, 9, 5, 6, 6, 11, 9, 1, 11, 8, 9, 3 ),1295( 3, 6, 11, 6, 0, 6, 3, 4, 0, 5, 6, 2, 0, 9, 5, 3 ),1296( 0, 11, 8, 6, 0, 5, 11, 0, 0, 1, 5, 3, 5, 6, 11, 3 ),1297( 6, 11, 3, 3, 6, 3, 5, 4, 5, 3, 1, 6 ),1298( 1, 2, 10, 7, 9, 5, 11, 1, 9, 11, 8, 6, 11, 5, 6, 6 ),1299( 0, 11, 3, 6, 0, 6, 11, 2, 0, 9, 6, 1, 5, 6, 9, 5, 1, 2, 10, 7 ),1300( 11, 8, 5, 1, 11, 5, 6, 6, 8, 0, 5, 1, 10, 5, 2, 5, 0, 2, 5, 1 ),1301( 6, 11, 3, 3, 6, 3, 5, 4, 2, 10, 3, 5, 10, 5, 3, 1 ),1302( 5, 8, 9, 6, 5, 2, 8, 0, 5, 6, 2, 3, 3, 8, 2, 5 ),1303( 9, 5, 6, 3, 9, 6, 0, 4, 0, 6, 2, 6 ),1304( 1, 5, 8, 1, 1, 8, 0, 6, 5, 6, 8, 1, 3, 8, 2, 5, 6, 2, 8, 1 ),1305( 1, 5, 6, 3, 2, 1, 6, 5 ),1306( 1, 3, 6, 1, 1, 6, 10, 6, 3, 8, 6, 1, 5, 6, 9, 5, 8, 9, 6, 1 ),1307( 10, 1, 0, 3, 10, 0, 6, 4, 9, 5, 0, 5, 5, 6, 0, 1 ),1308( 0, 3, 8, 7, 5, 6, 10, 7 ), ( 10, 5, 6, 7 ),1309( 11, 5, 10, 6, 7, 5, 11, 5 ), ( 11, 5, 10, 6, 11, 7, 5, 3, 8, 3, 0, 7 ),1310( 5, 11, 7, 6, 5, 10, 11, 3, 1, 9, 0, 7 ),1311( 10, 7, 5, 6, 10, 11, 7, 3, 9, 8, 1, 5, 8, 3, 1, 3 ),1312( 11, 1, 2, 6, 11, 7, 1, 1, 7, 5, 1, 3 ),1313( 0, 8, 3, 7, 1, 2, 7, 1, 1, 7, 5, 6, 7, 2, 11, 6 ),1314( 9, 7, 5, 6, 9, 2, 7, 0, 9, 0, 2, 3, 2, 11, 7, 3 ),1315( 7, 5, 2, 1, 7, 2, 11, 6, 5, 9, 2, 1, 3, 2, 8, 5, 9, 8, 2, 1 ),1316( 2, 5, 10, 6, 2, 3, 5, 1, 3, 7, 5, 3 ),1317( 8, 2, 0, 6, 8, 5, 2, 0, 8, 7, 5, 3, 10, 2, 5, 5 ),1318( 9, 0, 1, 7, 5, 10, 3, 1, 5, 3, 7, 6, 3, 10, 2, 6 ),1319( 9, 8, 2, 1, 9, 2, 1, 6, 8, 7, 2, 1, 10, 2, 5, 5, 7, 5, 2, 1 ),1320( 1, 3, 5, 5, 3, 7, 5, 3 ), ( 0, 8, 7, 3, 0, 7, 1, 4, 1, 7, 5, 6 ),1321( 9, 0, 3, 3, 9, 3, 5, 4, 5, 3, 7, 6 ), ( 9, 8, 7, 3, 5, 9, 7, 5 ),1322( 5, 8, 4, 6, 5, 10, 8, 1, 10, 11, 8, 3 ),1323( 5, 0, 4, 6, 5, 11, 0, 0, 5, 10, 11, 3, 11, 3, 0, 3 ),1324( 0, 1, 9, 7, 8, 4, 10, 1, 8, 10, 11, 6, 10, 4, 5, 6 ),1325( 10, 11, 4, 1, 10, 4, 5, 6, 11, 3, 4, 1, 9, 4, 1, 5, 3, 1, 4, 1 ),1326( 2, 5, 1, 6, 2, 8, 5, 0, 2, 11, 8, 3, 4, 5, 8, 5 ),1327( 0, 4, 11, 1, 0, 11, 3, 6, 4, 5, 11, 1, 2, 11, 1, 5, 5, 1, 11, 1 ),1328( 0, 2, 5, 1, 0, 5, 9, 6, 2, 11, 5, 1, 4, 5, 8, 5, 11, 8, 5, 1 ),1329( 9, 4, 5, 7, 2, 11, 3, 7 ),1330( 2, 5, 10, 6, 3, 5, 2, 4, 3, 4, 5, 2, 3, 8, 4, 3 ),1331( 5, 10, 2, 3, 5, 2, 4, 4, 4, 2, 0, 6 ),1332( 3, 10, 2, 6, 3, 5, 10, 2, 3, 8, 5, 1, 4, 5, 8, 5, 0, 1, 9, 7 ),1333( 5, 10, 2, 3, 5, 2, 4, 4, 1, 9, 2, 5, 9, 4, 2, 1 ),1334( 8, 4, 5, 3, 8, 5, 3, 4, 3, 5, 1, 6 ), ( 0, 4, 5, 3, 1, 0, 5, 5 ),1335( 8, 4, 5, 3, 8, 5, 3, 4, 9, 0, 5, 5, 0, 3, 5, 1 ), ( 9, 4, 5, 7 ),1336( 4, 11, 7, 6, 4, 9, 11, 1, 9, 10, 11, 3 ),1337( 0, 8, 3, 7, 4, 9, 7, 5, 9, 11, 7, 2, 9, 10, 11, 3 ),1338( 1, 10, 11, 3, 1, 11, 4, 0, 1, 4, 0, 6, 7, 4, 11, 5 ),1339( 3, 1, 4, 1, 3, 4, 8, 6, 1, 10, 4, 1, 7, 4, 11, 5, 10, 11, 4, 1 ),1340( 4, 11, 7, 6, 9, 11, 4, 4, 9, 2, 11, 2, 9, 1, 2, 3 ),1341( 9, 7, 4, 6, 9, 11, 7, 2, 9, 1, 11, 1, 2, 11, 1, 5, 0, 8, 3, 7 ),1342( 11, 7, 4, 3, 11, 4, 2, 4, 2, 4, 0, 6 ),1343( 11, 7, 4, 3, 11, 4, 2, 4, 8, 3, 4, 5, 3, 2, 4, 1 ),1344( 2, 9, 10, 6, 2, 7, 9, 0, 2, 3, 7, 3, 7, 4, 9, 3 ),1345( 9, 10, 7, 1, 9, 7, 4, 6, 10, 2, 7, 1, 8, 7, 0, 5, 2, 0, 7, 1 ),1346( 3, 7, 10, 1, 3, 10, 2, 6, 7, 4, 10, 1, 1, 10, 0, 5, 4, 0, 10, 1 ),1347( 1, 10, 2, 7, 8, 7, 4, 7 ), ( 4, 9, 1, 3, 4, 1, 7, 4, 7, 1, 3, 6 ),1348( 4, 9, 1, 3, 4, 1, 7, 4, 0, 8, 1, 5, 8, 7, 1, 1 ),1349( 4, 0, 3, 3, 7, 4, 3, 5 ), ( 4, 8, 7, 7 ),1350( 9, 10, 8, 5, 10, 11, 8, 3 ), ( 3, 0, 9, 3, 3, 9, 11, 4, 11, 9, 10, 6 ),1351( 0, 1, 10, 3, 0, 10, 8, 4, 8, 10, 11, 6 ),1352( 3, 1, 10, 3, 11, 3, 10, 5 ), ( 1, 2, 11, 3, 1, 11, 9, 4, 9, 11, 8, 6 ),1353( 3, 0, 9, 3, 3, 9, 11, 4, 1, 2, 9, 5, 2, 11, 9, 1 ),1354( 0, 2, 11, 3, 8, 0, 11, 5 ), ( 3, 2, 11, 7 ),1355( 2, 3, 8, 3, 2, 8, 10, 4, 10, 8, 9, 6 ), ( 9, 10, 2, 3, 0, 9, 2, 5 ),1356( 2, 3, 8, 3, 2, 8, 10, 4, 0, 1, 8, 5, 1, 10, 8, 1 ), ( 1, 10, 2, 7 ),1357( 1, 3, 8, 3, 9, 1, 8, 5 ), ( 0, 9, 1, 7 ), ( 0, 3, 8, 7 ), None )135813591360