Changeset 7906:115bdd54c299
- Timestamp:
- 01/02/08 18:45:26 (5 years ago)
- Branch:
- default
- Parents:
- 7887:2a12de9961da (diff), 7905:6645b83a4cd8 (diff)
Note: this is a merge changeset, the changes displayed below correspond to the merge itself.
Use the (diff) links above to see all the changes relative to each parent. - Location:
- sage
- Files:
-
- 11 edited
-
plot/plot3d/all.py (modified) (1 diff)
-
plot/plot3d/all.py (modified) (1 diff)
-
plot/plot3d/base.pyx (modified) (15 diffs)
-
plot/plot3d/base.pyx (modified) (1 diff)
-
plot/plot3d/index_face_set.pyx (modified) (6 diffs)
-
plot/plot3d/index_face_set.pyx (modified) (2 diffs)
-
plot/plot3d/parametric_surface.pyx (modified) (4 diffs)
-
plot/plot3d/parametric_surface.pyx (modified) (2 diffs)
-
plot/plot3d/plot3d.py (modified) (4 diffs)
-
server/notebook/guard.py (modified) (3 diffs)
-
server/notebook/guard.py (modified) (3 diffs)
Legend:
- Unmodified
- Added
- Removed
-
sage/plot/plot3d/all.py
r7886 r7906 1 1 from shapes import Box, ColorCube, Cone, Cylinder, LineSegment, Arrow, Sphere, Torus, Text as Text3D 2 2 from parametric_surface import ParametricSurface, MobiusStrip 3 from parametric_plot3d import parametric_plot3d 3 4 from plot3d import plot3d, axes as axes3d 4 5 from platonic import Tetrahedron, Cube, Octahedron, Dodecahedron, IndexFaceSet, Icosahedron 6 7 from shapes2 import line3d, frame3d 8 9 from texture import Texture -
sage/plot/plot3d/all.py
r7900 r7906 1 1 from shapes import Box, ColorCube, Cone, Cylinder, LineSegment, Arrow, Sphere, Torus, Text as Text3D 2 3 from parametric_surface import MobiusStrip 2 from parametric_surface import ParametricSurface, MobiusStrip 4 3 from parametric_plot3d import parametric_plot3d 5 6 4 from plot3d import plot3d, axes as axes3d 7 5 from platonic import Tetrahedron, Cube, Octahedron, Dodecahedron, IndexFaceSet, Icosahedron -
sage/plot/plot3d/base.pyx
r7887 r7906 75 75 return other 76 76 return Graphics3dGroup([self, other]) 77 78 def bounding_box(self): 79 raise NotImplementedError 77 80 78 81 def transform(self, **kwds): … … 130 133 resolution 400 400 131 134 132 133 135 camera 134 136 zoom 1.0 … … 136 138 antialiasing 1 137 139 raydepth 8 138 center 2. 0 1.0 0.75139 viewdir -2. 0 -1.0 -0.5140 center 2.3 2.4 2.0 141 viewdir -2.3 -2.4 -2.0 140 142 updir 0.0 0.0 1.0 141 143 end_camera … … 148 150 plane 149 151 center -2000 -1000 -500 150 normal -2.0 -1.0 -0.5152 normal 2.3 2.4 2.0 151 153 TEXTURE 152 154 AMBIENT 1.0 DIFFUSE 1.0 SPECULAR 1.0 OPACITY 1.0 … … 177 179 return "\n".join(flatten_list([self.obj_repr(render_params), ""])) 178 180 179 def export_jmol(self, filename='jmol_shape.jmol', force_reload=False, zoom=100, spin=False, background=(1,1,1), stereo=False): 181 def export_jmol(self, filename='jmol_shape.jmol', force_reload=False, 182 zoom=100, spin=False, background=(1,1,1), stereo=False, 183 perspective_depth = True, 184 orientation = (-764,-346,-545,76.39)): 185 # orientation chosen to look same as tachyon 180 186 render_params = self.default_render_params() 181 187 render_params.output_file = filename … … 191 197 if stereo is True: stereo = "redblue" 192 198 f.write('stereo %s\n' % stereo) 199 if orientation: 200 f.write('moveto 0 %s %s %s %s\n'%tuple(orientation)) 193 201 194 202 f.write('zoom %s\n'%zoom) 203 204 if perspective_depth: 205 f.write('set perspectivedepth ON\n') 206 else: 207 f.write('set perspectivedepth OFF\n') 195 208 196 209 # Put the rest of the object in … … 212 225 else: 213 226 return self.transform(T=T) 214 215 def show(self, viewer="jmol", filename="shape", verbosity=0, figsize=4, **kwds): 227 228 def _rescale_for_aspect_ratio_and_zoom(self, b, aspect_ratio, zoom): 229 if aspect_ratio is None: 230 return (b*zoom,b*zoom,b*zoom), (-b*zoom,-b*zoom,-b*zoom) 231 box = [b*w for w in aspect_ratio] 232 # Now take the maximum length in box and rescale to b. 233 s = b / max(box) 234 box_max = tuple([s*w*zoom for w in box]) 235 box_min = tuple([-w*zoom for w in box_max]) 236 return box_min, box_max 237 238 def _prepare_for_jmol(self, frame, axes, aspect_ratio, zoom): 239 box_min, box_max = self._rescale_for_aspect_ratio_and_zoom(6.0, aspect_ratio, zoom) 240 return self._transform_to_bounding_box(box_min, box_max, frame=frame, 241 axes=axes, thickness=1) 242 243 def _prepare_for_tachyon(self, frame, axes, aspect_ratio, zoom): 244 box_min, box_max = self._rescale_for_aspect_ratio_and_zoom(1.0, aspect_ratio, zoom) 245 A = self._transform_to_bounding_box(box_min, box_max, 246 frame=frame, axes=axes, thickness=0.5) 247 return A 248 249 def _transform_to_bounding_box(self, xyz_min, xyz_max, frame, axes, thickness): 250 a_min, a_max = self.bounding_box() 251 a_min = list(a_min); a_max = list(a_max) 252 for i in range(3): 253 if a_min[i] == a_max[i]: 254 a_min[i] = -1 255 a_max[i] = 1 256 257 # Rescale in each direction 258 scale = [(xyz_max[i] - xyz_min[i]) / (a_max[i] - a_min[i]) for i in range(3)] 259 X = self.scale(scale) 260 a_min = [scale[i]*a_min[i] for i in range(3)] 261 a_max = [scale[i]*a_max[i] for i in range(3)] 262 263 # Translate so lower left corner of original bounding box 264 # is in the right spot 265 T = [xyz_min[i] - a_min[i] for i in range(3)] 266 X = X.translate(T) 267 if frame: 268 from shapes2 import frame3d 269 F = frame3d(xyz_min, xyz_max, opacity=0.5, color=(0,0,0), thickness=thickness) 270 X += F 271 272 if axes: 273 # draw axes 274 from shapes import Arrow 275 A = (Arrow((min(0,a_min[0]),0, 0), (max(0,a_max[0]), 0,0), 276 thickness, color="blue"), 277 Arrow((0,min(0,a_min[1]), 0), (0, max(0,a_max[1]), 0), 278 thickness, color="blue"), 279 Arrow((0, 0, min(0,a_min[2])), (0, 0, max(0,a_max[2])), 280 thickness, color="blue")) 281 X += sum(A).translate([-z for z in T]) 282 283 return X 284 285 def show(self, viewer="jmol", filename=None, verbosity=0, figsize=5, 286 aspect_ratio = None, zoom=1, 287 frame=True, axes = False, **kwds): 216 288 """ 217 289 INPUT: … … 220 292 'tachyon': an embedded ray tracer 221 293 'java3d': a popup OpenGL 3d java applet 222 filename -- string (default: 'shape'); file to save the image to294 filename -- string (default: a temp file); file to save the image to 223 295 verbosity -- display information about rendering the figure 224 figsize -- (default: 4); x or pair [x,y] for numbers, e.g., [4,4]; controls225 the size of the output figure. E.g., with jmolthe number of296 figsize -- (default: 5); x or pair [x,y] for numbers, e.g., [5,5]; controls 297 the size of the output figure. E.g., with Tachyon the number of 226 298 pixels in each direction is 100 times figsize[0]. 299 This is ignored for the jmol embedded renderer. 227 300 **kwds -- other options, which make sense for particular rendering engines 228 301 """ 302 import sage.misc.misc 303 if filename is None: 304 filename = sage.misc.misc.tmp_filename() 229 305 if not isinstance(figsize, (list,tuple)): 230 306 figsize = [figsize, figsize] 231 307 from sage.plot.plot import EMBEDDED_MODE, DOCTEST_MODE 232 import sage.misc.misc233 308 ext = None 309 310 # Tachyon resolution options 234 311 if DOCTEST_MODE: 235 312 opts = '-res 10 10' 236 313 filename = sage.misc.misc.SAGE_TMP + "/tmp" 237 314 elif EMBEDDED_MODE: 238 opts = ' '315 opts = '-res %s %s'%(figsize[0]*100, figsize[1]*100) 239 316 filename = sage.misc.misc.graphics_filename()[:-4] 240 317 else: 241 opts = ' '318 opts = '-res %s %s'%(figsize[0]*100, figsize[1]*100) 242 319 243 320 if DOCTEST_MODE or viewer=='tachyon' or (viewer=='java3d' and EMBEDDED_MODE): 244 tachyon_rt(self.tachyon(**kwds), filename+".png", verbosity, True, opts) 321 T = self._prepare_for_tachyon(frame, axes, aspect_ratio, zoom) 322 tachyon_rt(T.tachyon(**kwds), filename+".png", verbosity, True, opts) 245 323 ext = "png" 246 324 import sage.misc.viewer 247 325 viewer_app = sage.misc.viewer.browser() 326 248 327 if DOCTEST_MODE or viewer=='java3d': 249 328 f = open(filename+".obj", "w") … … 258 337 259 338 if DOCTEST_MODE or viewer=='jmol': 260 # Encode the desired applet size in the end of the filename: 339 # Temporary hack: encode the desired applet size in the end of the filename: 340 # (This will be removed once we have dynamic resizing of applets in the browser.) 261 341 base, ext = os.path.splitext(filename) 262 filename = '%s-size%s%s'%(base, figsize[0]*100, ext) 263 self.export_jmol(filename + ".jmol", force_reload=EMBEDDED_MODE, **kwds) 342 fg = figsize[0] 343 #if fg >= 2: 344 # fg = 2 345 filename = '%s-size%s%s'%(base, fg*100, ext) 346 347 T = self._prepare_for_jmol(frame, axes, aspect_ratio, zoom) 348 T.export_jmol(filename + ".jmol", force_reload=EMBEDDED_MODE, **kwds) 264 349 viewer_app = sage.misc.misc.SAGE_LOCAL + "/java/jmol/jmol" 265 350 ext = "jmol" … … 267 352 if ext is None: 268 353 raise ValueError, "Unknown 3d plot type: %s" % viewer 354 269 355 if not DOCTEST_MODE and not EMBEDDED_MODE: 270 356 if verbosity: … … 277 363 def __init__(self, all=[]): 278 364 self.all = all 365 366 def bounding_box(self): 367 # Box that contains the bounding boxes of 368 # all the objects that make up self. 369 v = [obj.bounding_box() for obj in self.all] 370 return min3([a[0] for a in v]), max3([a[1] for a in v]) 279 371 280 372 def transform(self, **kwds): … … 323 415 if T is not None: 324 416 self.T = T 417 418 def bounding_box(self): 419 try: 420 return self._bounding_box 421 except AttributeError: 422 pass 423 # Get the box before transformation 424 a = Graphics3dGroup.bounding_box(self) 425 # The corners of the box 426 import sage.misc.mrange 427 corners = [] 428 for f in sage.misc.mrange.cartesian_product_iterator([[0,1]]*3): 429 corners.append([a[f[i]][i] for i in range(3)]) 430 # Transform the corners of the box 431 T = self.get_transformation() 432 w = [T.transform_point(p) for p in corners] 433 # Figure out what the new bounding box is 434 a_min = [min([z[i] for z in w]) for i in range(3)] 435 a_max = [max([z[i] for z in w]) for i in range(3)] 436 self._bounding_box = a_min, a_max 437 return self._bounding_box 325 438 326 439 def transform(self, **kwds): … … 392 505 393 506 394 cdef class Prim ativeObject(Graphics3d):507 cdef class PrimitiveObject(Graphics3d): 395 508 def __init__(self, **kwds): 396 509 try: … … 402 515 self.texture = kwds['color'] 403 516 if not is_Texture(self.texture): 404 self.texture = Texture(self.texture) 517 if kwds.has_key('opacity'): 518 self.texture = Texture(self.texture, opacity=kwds['opacity']) 519 else: 520 self.texture = Texture(self.texture) 405 521 except KeyError: 406 522 self.texture = default_texture 407 523 408 def set_texture(self, texture ):524 def set_texture(self, texture, **kwds): 409 525 if not is_Texture(texture): 410 texture = Texture(texture )526 texture = Texture(texture, **kwds) 411 527 self.texture = texture 412 528 … … 517 633 i += 1 518 634 return flat 635 636 637 def min3(v): 638 """ 639 Return the componentwise minimum of a list of 3-tuples. 640 641 EXAMPLES: 642 sage: from sage.plot.plot3d.base import min3, max3 643 sage: min3([(-1,2,5), (-3, 4, 2)]) 644 (-3, 2, 2) 645 """ 646 return tuple([min([a[i] for a in v]) for i in range(3)]) 647 648 def max3(v): 649 """ 650 Return the componentwise maximum of a list of 3-tuples. 651 652 EXAMPLES: 653 sage: from sage.plot.plot3d.base import min3, max3 654 sage: max3([(-1,2,5), (-3, 4, 2)]) 655 (-1, 4, 5) 656 """ 657 return tuple([max([a[i] for a in v]) for i in range(3)]) -
sage/plot/plot3d/base.pyx
r7905 r7906 581 581 self.transform = None 582 582 self.ds = 1 583 self.crease_threshold = .8 583 584 self.__dict__.update(kwds) 584 585 -
sage/plot/plot3d/index_face_set.pyx
r7886 r7906 1 """ nodoctest1 """ 2 2 Graphics3D object that consists of a list of polygons, also used for 3 3 triangulations of other objects. … … 144 144 145 145 146 cdef class IndexFaceSet(Prim ativeObject):146 cdef class IndexFaceSet(PrimitiveObject): 147 147 148 148 """ … … 178 178 """ 179 179 180 def __new__(self, faces, point_list=None, enclos de=False, **kwds):180 def __new__(self, faces, point_list=None, enclosed=False, **kwds): 181 181 self.vs = <point_c *>NULL 182 182 self.face_indices = <int *>NULL … … 185 185 186 186 def __init__(self, faces, point_list=None, enclosed=False, **kwds): 187 Prim ativeObject.__init__(self, **kwds)187 PrimitiveObject.__init__(self, **kwds) 188 188 189 189 self.enclosed = enclosed … … 464 464 465 465 def bounding_box(self): 466 if self.vcount == 0: 467 return ((0,0,0),(0,0,0)) 468 466 469 cdef Py_ssize_t i 467 470 cdef point_c low = self.vs[0], high = self.vs[0] … … 647 650 f.write('\n') 648 651 f.close() 649 if render_params.force_reload: 650 filename += "?%s" % randint(1,1000000) 651 return ['pmesh %s "%s"\n%s' % (name, filename, self.texture.jmol_str("pmesh"))] 652 #if render_params.force_reload: 653 # filename += "?%s" % randint(1,1000000) 654 655 s = 'pmesh %s "%s"\n%s' % (name, filename, self.texture.jmol_str("pmesh")) 656 657 # If we wanted to turn on display of the mesh lines or dots 658 # we would uncomment thse. This should be determined by 659 # render_params, probably. 660 #s += '\npmesh %s mesh\n'%name 661 #s += '\npmesh %s dots\n'%name 662 return [s] 652 663 653 664 def dual(self, **kwds): -
sage/plot/plot3d/index_face_set.pyx
r7900 r7906 271 271 self.vcount = ix 272 272 sage_free(point_map) 273 274 def _seperate_creases(self, threshold): 275 """ 276 Some rendering engines Gouraud shading, which is great for smooth 277 surfaces but looks bad if one actually has a polyhedron. 278 279 INPUT: 280 threshold -- the minimum cosine of the angle between adjacent faces 281 a higher threshold seperates more, all faces if >= 1, 282 no faces if <= -1 283 """ 284 cdef Py_ssize_t i, j, k 285 cdef face_c *face 286 cdef int v, count, total = 0 287 cdef int* point_counts = <int *>sage_malloc(sizeof(int) * (self.vcount * 2 + 1)) 288 # For each vertex, get number of faces 289 if point_counts == NULL: 290 raise MemoryError, "Out of memory in _seperate_creases for %s" % type(self) 291 cdef int* running_point_counts = &point_counts[self.vcount] 292 memset(point_counts, 0, sizeof(int) * self.vcount) 293 for i from 0 <= i < self.fcount: 294 face = &self._faces[i] 295 total += face.n 296 for j from 0 <= j < face.n: 297 point_counts[face.vertices[j]] += 1 298 # Running used as index into face list 299 cdef int running = 0 300 cdef int max = 0 301 for i from 0 <= i < self.vcount: 302 running_point_counts[i] = running 303 running += point_counts[i] 304 if point_counts[i] > max: 305 max = point_counts[i] 306 running_point_counts[self.vcount] = running 307 # Create an array, indexed by running_point_counts[v], to the list of faces containing that vertex. 308 cdef face_c** point_faces = <face_c **>sage_malloc(sizeof(face_c*) * total) 309 if point_faces == NULL: 310 sage_free(point_counts) 311 raise MemoryError, "Out of memory in _seperate_creases for %s" % type(self) 312 _sig_on 313 memset(point_counts, 0, sizeof(int) * self.vcount) 314 for i from 0 <= i < self.fcount: 315 face = &self._faces[i] 316 for j from 0 <= j < face.n: 317 v = face.vertices[j] 318 point_faces[running_point_counts[v]+point_counts[v]] = face 319 point_counts[v] += 1 320 # Now, for each vertex, see if all faces are close enough, 321 # or if it is a crease. 322 cdef face_c** faces 323 cdef int start = 0 324 cdef bint any 325 # We compare against face 0, and if it's not flat enough we push it to the end. 326 # Then we come around again to compare everything that was put at the end, possibly 327 # pushing stuff to the end again (until no further changes are needed). 328 while start < self.vcount: 329 ix = self.vcount 330 # Find creases 331 for i from 0 <= i < self.vcount - start: 332 faces = &point_faces[running_point_counts[i]] 333 any = 0 334 for j from point_counts[i] > j >= 1: 335 if cos_face_angle(faces[0][0], faces[j][0], self.vs) < threshold: 336 any = 1 337 face = faces[j] 338 point_counts[i] -= 1 339 if j != point_counts[i]: 340 faces[j] = faces[point_counts[i]] # swap 341 faces[point_counts[i]] = face 342 if any: 343 ix += 1 344 # Reallocate room for vertices at end 345 if ix > self.vcount: 346 self.vs = <point_c *>sage_realloc(self.vs, sizeof(point_c) * ix) 347 if self.vs == NULL: 348 sage_free(point_counts) 349 sage_free(point_faces) 350 self.vcount = self.fcount = self.icount = 0 # so we don't get segfaults on bad points 351 _sig_off 352 raise MemoryError, "Out of memory in _seperate_creases for %s, CORRUPTED" % type(self) 353 ix = self.vcount 354 running = 0 355 for i from 0 <= i < self.vcount - start: 356 if point_counts[i] != running_point_counts[i+1] - running_point_counts[i]: 357 # We have a new vertex 358 self.vs[ix] = self.vs[i+start] 359 # Update the point_counts and point_faces arrays for the next time around. 360 count = running_point_counts[i+1] - running_point_counts[i] - point_counts[i] 361 faces = &point_faces[running] 362 for j from 0 <= j < count: 363 faces[j] = point_faces[running_point_counts[i] + point_counts[i] + j] 364 face = faces[j] 365 for k from 0 <= k < face.n: 366 if face.vertices[k] == i + start: 367 face.vertices[k] = ix 368 point_counts[ix-self.vcount] = count 369 running_point_counts[ix-self.vcount] = running 370 running += count 371 ix += 1 372 running_point_counts[ix-self.vcount] = running 373 start = self.vcount 374 self.vcount = ix 375 376 sage_free(point_counts) 377 sage_free(point_faces) 378 _sig_off 379 380 273 381 274 382 def _mem_stats(self): … … 506 614 cdef Py_ssize_t i 507 615 cdef point_c res 508 616 617 self._seperate_creases(render_params.crease_threshold) 618 509 619 _sig_on 510 620 if transform is None: -
sage/plot/plot3d/parametric_surface.pyx
r7886 r7906 79 79 sage: len(S.face_list()) 80 80 2214 81 82 The Hessenberg surface: 83 sage: def f(u,v): 84 ... a = 1 85 ... from math import cos, sin, sinh, cosh 86 ... x = cos(a)*(cos(u)*sinh(v)-cos(3*u)*sinh(3*v)/3)+ \ 87 ... sin(a)*(sin(u)*cosh(v)-sin(3*u)*cosh(3*v)/3) 88 ... y = cos(a)*(sin(u)*sinh(v)+sin(3*u)*sinh(3*v)/3)+\ 89 ... sin(a)*(-cos(u)*cosh(v)-cos(3*u)*cosh(3*v)/3) 90 ... z = cos(a)*cos(2*u)*cosh(2*v)+sin(a)*sin(2*u)*sinh(2*v) 91 ... return (x,y,z) 92 sage: v = srange(float(0),float((3/2)*pi),float(0.1)) 93 sage: S = ParametricSurface(f, (srange(float(0),float(pi),float(0.1)), \ 94 ... srange(float(-1),float(1),float(0.1))), color="blue") 95 sage: show(S) 81 96 """ 82 97 … … 127 142 self.triangulate() 128 143 return IndexFaceSet.dual(self) 144 145 def bounding_box(self): 146 # We must triangulate before computing the bounding box; otherwise 147 # we'll get an empty bounding box, as the bounding box is computed 148 # using the triangulation, and before triangulating the triangulation 149 # is empty. 150 self.triangulate() 151 return IndexFaceSet.bounding_box(self) 129 152 130 153 def triangulate(self, render_params=None): … … 149 172 # Already triangulated at on this grid. 150 173 return 151 174 152 175 cdef Py_ssize_t i, j 153 176 cdef Py_ssize_t n = len(urange) - 1 … … 163 186 self.eval_c(&self.vs[ix], u, v) 164 187 ix += 1 165 except: 188 except: # TODO -- this would catch control-C,etc. -- FIX THIS TO CATCH WHAT IS RAISED!!!! 166 189 _sig_off 167 190 self.fcount = self.vcount = 0 -
sage/plot/plot3d/parametric_surface.pyx
r7900 r7906 55 55 from math import cos, sin 56 56 from sage.rings.all import RDF 57 58 from base import RenderParams 57 59 58 60 … … 98 100 self.render_grid = domain 99 101 IndexFaceSet.__init__(self, [], [], **kwds) 102 103 def default_render_params(self): 104 return RenderParams(ds=.075, crease_threshold=.35) 100 105 101 106 def tachyon_repr(self, render_params): -
sage/plot/plot3d/plot3d.py
r7886 r7906 60 60 61 61 62 def plot3d(f,(xmin,xmax),(ymin,ymax),texture=None, grad_f=None,63 max_bend=. 7,max_depth=6,initial_depth=4, num_colors=None):62 def plot3d(f,(xmin,xmax),(ymin,ymax),texture=None, opacity=1, grad_f=None, 63 max_bend=.5, max_depth=5, initial_depth=4, num_colors=None): 64 64 """ 65 65 EXAMPLES: … … 67 67 68 68 """ 69 if initial_depth >= max_depth: 70 max_depth = initial_depth 71 xmin = float(xmin) 72 xmax = float(xmax) 73 ymin = float(ymin) 74 ymax = float(ymax) 75 69 76 # Check if f has a fast float evaluation 70 77 #try: … … 73 80 # # Nope -- no prob. 74 81 # pass 82 if texture is None: 83 texture = rainbow(128, 'rgbtuple') 75 84 76 if texture is None:77 texture = rainbow(100, 'rgbtuple')78 85 factory = TrivialTriangleFactory() 79 86 plot = TrianglePlot(factory, f, (xmin, xmax), (ymin, ymax), g = grad_f, 80 min_depth=initial_depth, max_depth=max_depth, max_bend=max_bend, num_colors = num_colors) 87 min_depth=initial_depth, max_depth=max_depth, 88 max_bend=max_bend, num_colors = num_colors) 81 89 82 90 P = IndexFaceSet(plot._objects) … … 92 100 min_z = bounds[0][2] 93 101 max_z = bounds[1][2] 94 span = len(texture) / (max_z - min_z) 102 if max_z == min_z: 103 span = 0 104 else: 105 span = len(texture) / (max_z - min_z) # max to avoid dividing by 0 95 106 parts = P.partition(lambda x,y,z: int((z-min_z)*span)) 96 107 all = [] 97 108 for k, G in parts.iteritems(): 98 G.set_texture(texture[k] )109 G.set_texture(texture[k], opacity=opacity) 99 110 all.append(G) 100 111 return Graphics3dGroup(all) -
sage/server/notebook/guard.py
r5208 r7906 150 150 authentication status of a given user. 151 151 """ 152 #log.msg("request: %s, segments: %s" % (str(request.args), segments))152 log.msg("request: %s, segments: %s" % (str(request.args), segments)) 153 153 #see if the user already has a session going 154 154 if segments and segments[0] == "login": … … 160 160 l.addCallback(lambda _: self.requestPasswordAuthentication(request, segments)) 161 161 return l 162 162 163 163 session = self.getSession(request) 164 164 if session is None: 165 # log.msg("unknown session")165 # log.msg("unknown session") 166 166 return self.requestAnonymousAuthentication(request, segments) 167 167 else: … … 170 170 return self.logout(session, request, segments) 171 171 else: 172 #log.msg("session found ... locateResource")172 log.msg("session found ... locateResource") 173 173 creds = session.get_authCreds() 174 174 return self.locateResource(request, segments, session, creds) -
sage/server/notebook/guard.py
r5208 r7906 150 150 authentication status of a given user. 151 151 """ 152 #log.msg("request: %s, segments: %s" % (str(request.args), segments))152 log.msg("request: %s, segments: %s" % (str(request.args), segments)) 153 153 #see if the user already has a session going 154 154 if segments and segments[0] == "login": … … 160 160 l.addCallback(lambda _: self.requestPasswordAuthentication(request, segments)) 161 161 return l 162 162 163 163 session = self.getSession(request) 164 164 if session is None: 165 # log.msg("unknown session")165 # log.msg("unknown session") 166 166 return self.requestAnonymousAuthentication(request, segments) 167 167 else: … … 170 170 return self.logout(session, request, segments) 171 171 else: 172 #log.msg("session found ... locateResource")172 log.msg("session found ... locateResource") 173 173 creds = session.get_authCreds() 174 174 return self.locateResource(request, segments, session, creds)
Note: See TracChangeset
for help on using the changeset viewer.
