def test_edge_curves_elasticity(self): # create an arrow-like 2D geometry with the pointy end at (-1,1) towards up and left # rebuild to avoid rational representations c1 = CurveFactory.circle_segment(pi / 2).rebuild(3,11) c2 = Curve(BSplineBasis(2, [0, 0, 1, 2, 2]), [[0, 1], [-1, 1], [-1, 0]]) c3 = CurveFactory.circle_segment(pi / 2).rebuild(3,11) c3.rotate(pi) c4 = Curve(BSplineBasis(2), [[0, -1], [1, 0]]).rebuild(3,10) c4 = c4.rebuild(4,11) surf = SurfaceFactory.edge_curves([c1, c2, c3, c4], type='elasticity') # check right dimensions of output self.assertEqual(surf.shape[0], 11) # 11 controlpoints in the circle segment self.assertEqual(surf.shape[1], 13) # 11 controlpoints in c4, +2 for C0-knot in c1 self.assertEqual(surf.order(0), 3) self.assertEqual(surf.order(1), 4) # check that c1 edge conforms to surface edge u = np.linspace(0,1,7) c1.reparam() pts_surf = surf(u,0.0) pts_c1 = c1(u) for (xs,xc) in zip(pts_surf[:,0,:], pts_c1): self.assertTrue(np.allclose(xs, xc)) # check that c2 edge conforms to surface edge v = np.linspace(0,1,7) c2.reparam() pts_surf = surf(1.0,v) pts_c2 = c2(v) for (xs,xc) in zip(pts_surf[:,0,:], pts_c2): self.assertTrue(np.allclose(xs, xc))
def test_edge_curves(self): # create an arrow-like 2D geometry with the pointy end at (-1,1) towards up and left # mixes rational and non-rational curves with different parametrization spaces c1 = cf.circle_segment(pi / 2) c2 = Curve(BSplineBasis(2, [0, 0, 1, 2, 2]), [[0, 1], [-1, 1], [-1, 0]]) c3 = cf.circle_segment(pi / 2) c3.rotate(pi) c4 = Curve(BSplineBasis(2), [[0, -1], [1, 0]]) surf = sf.edge_curves(c1, c2, c3, c4) # srf spits out parametric space (0,1)^2, so we sync these up to input curves c3.reverse() c4.reverse() c1.reparam() c2.reparam() c3.reparam() c4.reparam() for u in np.linspace(0, 1, 7): self.assertAlmostEqual(surf(u, 0)[0], c1(u)[0]) # x-coord, bottom crv self.assertAlmostEqual(surf(u, 0)[1], c1(u)[1]) # y-coord, bottom crv for u in np.linspace(0, 1, 7): self.assertAlmostEqual(surf(u, 1)[0], c3(u)[0]) # x-coord, top crv self.assertAlmostEqual(surf(u, 1)[1], c3(u)[1]) # y-coord, top crv for v in np.linspace(0, 1, 7): self.assertAlmostEqual(surf(0, v)[0], c4(v)[0]) # x-coord, left crv self.assertAlmostEqual(surf(0, v)[1], c4(v)[1]) # y-coord, left crv for v in np.linspace(0, 1, 7): self.assertAlmostEqual(surf(1, v)[0], c2(v)[0]) # x-coord, right crv self.assertAlmostEqual(surf(1, v)[1], c2(v)[1]) # y-coord, right crv # add a case where opposing sites have mis-matching rationality crvs = Surface().edges() # returned in order umin, umax, vmin, vmax crvs[0].force_rational() crvs[1].reverse() crvs[2].reverse() # input curves should be clockwise oriented closed loop srf = sf.edge_curves(crvs[0], crvs[3], crvs[1], crvs[2]) crvs[1].reverse() u = np.linspace(0, 1, 7) self.assertTrue(np.allclose(srf(u, 0).reshape((7, 2)), crvs[0](u))) self.assertTrue(np.allclose(srf(u, 1).reshape((7, 2)), crvs[1](u))) # test self-organizing curve ordering when they are not sequential srf = sf.edge_curves(crvs[0], crvs[2].reverse(), crvs[3], crvs[1]) u = np.linspace(0, 1, 7) self.assertTrue(np.allclose(srf(u, 0).reshape((7, 2)), crvs[0](u))) self.assertTrue(np.allclose(srf(u, 1).reshape((7, 2)), crvs[1](u))) # test error handling with self.assertRaises(ValueError): srf = sf.edge_curves(crvs + (Curve(), )) # 5 input curves
def teapot(): """ Generate the Utah teapot as 32 cubic bezier patches. This teapot has a rim, but no bottom. It is also self-intersecting making it unsuitable for perfect-match multipatch modeling. The data is picked from http://www.holmes3d.net/graphics/teapot/ :return: The utah teapot :rtype: List of Surface """ path = join(dirname(realpath(__file__)), 'templates', 'teapot.bpt') with open(path) as f: results = [] numb_patches = int(f.readline()) for i in range(numb_patches): p = np.fromstring(f.readline(), dtype=np.uint8, count=2, sep=' ') basis1 = BSplineBasis(p[0] + 1) basis2 = BSplineBasis(p[1] + 1) ncp = basis1.num_functions() * basis2.num_functions() cp = [ np.fromstring(f.readline(), dtype=np.float, count=3, sep=' ') for j in range(ncp) ] results.append(Surface(basis1, basis2, cp)) return results
def test_integrate(self): # create the linear functions x(t) = [1-t, t] on t=[0,1] b = BSplineBasis() self.assertAlmostEqual(b.integrate(0,1)[0], 0.5) self.assertAlmostEqual(b.integrate(0,1)[1], 0.5) self.assertAlmostEqual(b.integrate(.25,.5)[0], 5.0/32) self.assertAlmostEqual(b.integrate(.25,.5)[1], 3.0/32) # create the quadratic functions x(t) = [(1-t)^2, 2t(1-t), t^2] on t=[0,1] b = BSplineBasis(3) self.assertAlmostEqual(b.integrate(0,1)[0], 1.0/3) self.assertAlmostEqual(b.integrate(0,1)[1], 1.0/3) self.assertAlmostEqual(b.integrate(0,1)[2], 1.0/3) self.assertAlmostEqual(b.integrate(.25,.5)[0], 19.0/192) self.assertAlmostEqual(b.integrate(.25,.5)[1], 11.0/96) self.assertAlmostEqual(b.integrate(.25,.5)[2], 7.0/192) # create periodic quadratic functions on [0,3]. This is 3 functions, which are all # translated versions of the one below: # | 1/2 t^2 t=[0,1] # N[3] = { -t^2 + 3t - 3/2 t=[1,2] # | 1/2 (3-t)^2 t=[2,3] b = BSplineBasis(3, [-2,-1,0,1,2,3,4,5], periodic=1) self.assertEqual(len(b.integrate(0,3)), 3) # returns 3 functions self.assertAlmostEqual(b.integrate(0,3)[0], 1) self.assertAlmostEqual(b.integrate(0,3)[1], 1) self.assertAlmostEqual(b.integrate(0,3)[2], 1) self.assertAlmostEqual(b.integrate(0,1)[0], 1.0/6) self.assertAlmostEqual(b.integrate(0,1)[1], 2.0/3) self.assertAlmostEqual(b.integrate(0,1)[2], 1.0/6) self.assertAlmostEqual(b.integrate(0,2)[0], 2.0/6) self.assertAlmostEqual(b.integrate(0,2)[1], 5.0/6) self.assertAlmostEqual(b.integrate(0,2)[2], 5.0/6)
def test_surface_loft(self): crv1 = Curve(BSplineBasis(3, range(11), 1), [[1,-1], [1,0], [1,1], [-1,1], [-1,0], [-1,-1]]) crv2 = cf.circle(2) + (0,0,1) crv3 = Curve(BSplineBasis(4, range(11), 2), [[1,-1,2], [1,1,2], [-1,1,2], [-1,-1,2]]) crv4 = cf.circle(2) + (0,0,3) surf = sf.loft(crv1, crv2, crv3, crv4) crv1.set_dimension(3) # for convenience when evaluating t = np.linspace( 0, 1, 13) u = np.linspace(crv1.start(0), crv1.end(0), 13) pt = crv1(u) pt2 = surf(t,0).reshape(13,3) self.assertAlmostEqual(np.linalg.norm(pt-pt2), 0.0) u = np.linspace(crv2.start(0), crv2.end(0), 13) pt = crv2(u) pt2 = surf(t,1).reshape(13,3) self.assertAlmostEqual(np.linalg.norm(pt-pt2), 0.0) u = np.linspace(crv3.start(0), crv3.end(0), 13) pt = crv3(u) pt2 = surf(t,2).reshape(13,3) self.assertAlmostEqual(np.linalg.norm(pt-pt2), 0.0) u = np.linspace(crv4.start(0), crv4.end(0), 13) pt = crv4(u) pt2 = surf(t,3).reshape(13,3) self.assertAlmostEqual(np.linalg.norm(pt-pt2), 0.0)
def test_start_end(self): b = BSplineBasis(4, [-2, -1, 0, 0, 1, 2, 3, 3, 4, 5], periodic=1) self.assertAlmostEqual(b.start(), 0) self.assertAlmostEqual(b.end(), 3) b = BSplineBasis(3, [0, 1, 1, 2, 3, 4, 4, 6]) self.assertAlmostEqual(b.start(), 1) self.assertAlmostEqual(b.end(), 4)
def test_insert_knot(self): # more or less random 2D surface with p=[3,2] and n=[4,3] controlpoints = [[0, 0], [-1, 1], [0, 2], [1, -1], [1, 0], [1, 1], [2, 1], [2, 2], [2, 3], [3, 0], [4, 1], [3, 2]] basis1 = BSplineBasis(4, [0, 0, 0, 0, 2, 2, 2, 2]) basis2 = BSplineBasis(3, [0, 0, 0, 1, 1, 1]) surf = Surface(basis1, basis2, controlpoints) # pick some evaluation point (could be anything) evaluation_point1 = surf(0.23, 0.37) surf.insert_knot(.22, 0) surf.insert_knot(.5, 0) surf.insert_knot(.7, 0) surf.insert_knot(.1, 1) surf.insert_knot(1.0 / 3, 1) knot1, knot2 = surf.knots(with_multiplicities=True) self.assertEqual(len(knot1), 11) # 8 to start with, 3 new ones self.assertEqual(len(knot2), 8) # 6 to start with, 2 new ones evaluation_point2 = surf(0.23, 0.37) # evaluation before and after InsertKnot should remain unchanged self.assertAlmostEqual(evaluation_point1[0], evaluation_point2[0]) self.assertAlmostEqual(evaluation_point1[1], evaluation_point2[1]) # test a rational 2D surface controlpoints = [[0, 0, 1], [-1, 1, .96], [0, 2, 1], [1, -1, 1], [1, 0, .8], [1, 1, 1], [2, 1, .89], [2, 2, .9], [2, 3, 1], [3, 0, 1], [4, 1, 1], [3, 2, 1]] basis1 = BSplineBasis(3, [0, 0, 0, .4, 1, 1, 1]) basis2 = BSplineBasis(3, [0, 0, 0, 1, 1, 1]) surf = Surface(basis1, basis2, controlpoints, True) evaluation_point1 = surf(0.23, 0.37) surf.insert_knot(.22, 0) surf.insert_knot(.5, 0) surf.insert_knot(.7, 0) surf.insert_knot(.1, 1) surf.insert_knot(1.0 / 3, 1) knot1, knot2 = surf.knots(with_multiplicities=True) self.assertEqual(len(knot1), 10) # 7 to start with, 3 new ones self.assertEqual(len(knot2), 8) # 6 to start with, 2 new ones evaluation_point2 = surf(0.23, 0.37) # evaluation before and after InsertKnot should remain unchanged self.assertAlmostEqual(evaluation_point1[0], evaluation_point2[0]) self.assertAlmostEqual(evaluation_point1[1], evaluation_point2[1]) # test errors and exceptions with self.assertRaises(TypeError): surf.insert_knot(1, 2, 3) # too many arguments with self.assertRaises(ValueError): surf.insert_knot("tree-fiddy", .5) # wrong argument type with self.assertRaises(ValueError): surf.insert_knot(0, -0.2) # Outside-domain error with self.assertRaises(ValueError): surf.insert_knot(1, 1.4) # Outside-domain error
def read_surface(self, nsrf): knotsu = [0] for i in nsrf.KnotsU: knotsu.append(i) knotsu.append(knotsu[len(knotsu) - 1]) knotsu[0] = knotsu[1] knotsv = [0] for i in nsrf.KnotsV: knotsv.append(i) knotsv.append(knotsv[len(knotsv) - 1]) knotsv[0] = knotsv[1] basisu = BSplineBasis(nsrf.OrderU, knotsu, -1) basisv = BSplineBasis(nsrf.OrderV, knotsv, -1) cpts = [] cpts = np.ndarray( (nsrf.Points.CountU * nsrf.Points.CountV, 3 + nsrf.IsRational)) for v in range(0, nsrf.Points.CountV): for u in range(0, nsrf.Points.CountU): cpts[u + v * nsrf.Points.CountU, 0] = nsrf.Points[u, v].X cpts[u + v * nsrf.Points.CountU, 1] = nsrf.Points[u, v].Y cpts[u + v * nsrf.Points.CountU, 2] = nsrf.Points[u, v].Z if nsrf.IsRational: cpts[u + v * nsrf.Points.CountU, 3] = nsrf.Points[u, v].W return Surface(basisu, basisv, cpts, nsrf.IsRational)
def test_make_periodic(self): my_cps = np.array([[0, -1], [1, 0], [0, 1], [-1, 0], [0, -1]], dtype=float) crv = Curve(BSplineBasis(2, [0, 0, 1, 2, 3, 4, 4]), my_cps, rational=False) crv = crv.make_periodic(0) cps = [[0, -1], [1, 0], [0, 1], [-1, 0]] self.assertAlmostEqual(np.linalg.norm(crv.controlpoints - cps), 0.0) crv = Curve(BSplineBasis(2, [0, 0, 1, 2, 3, 4, 4]), my_cps, rational=False) crv = crv.make_periodic(0) cps = [[0, -1], [1, 0], [0, 1], [-1, 0]] self.assertAlmostEqual(np.linalg.norm(crv.controlpoints - cps), 0.0) crv = Curve(BSplineBasis(3, [0, 0, 0, 1, 2, 3, 3, 3]), my_cps, rational=False) crv = crv.make_periodic(0) cps = [[0, -1], [1, 0], [0, 1], [-1, 0]] self.assertAlmostEqual(np.linalg.norm(crv.controlpoints - cps), 0.0) crv = Curve(BSplineBasis(3, [0, 0, 0, 1, 2, 3, 3, 3]), my_cps, rational=False) crv = crv.make_periodic(1) cps = [[-1, 0], [1, 0], [0, 1]] self.assertAlmostEqual(np.linalg.norm(crv.controlpoints - cps), 0.0)
def test_make_identical(self): basis1 = BSplineBasis(4, [-1, -1, 0, 0, 1, 2, 9, 9, 10, 10, 11, 12], periodic=1) basis2 = BSplineBasis(2, [0, 0, .25, 1, 1]) cp = [[0, 0, 0, 1], [1, 0, 0, 1.1], [2, 0, 0, 1], [0, 1, 0, .7], [1, 1, 0, .8], [2, 1, 0, 1], [0, 2, 0, 1], [1, 2, 0, 1.2], [2, 2, 0, 1]] surf1 = Surface(BSplineBasis(3), BSplineBasis(3), cp, True) # rational 3D surf2 = Surface(basis1, basis2) # periodic 2D surf1.insert_knot([0.25, .5, .75], direction='v') Surface.make_splines_identical(surf1, surf2) for s in (surf1, surf2): self.assertEqual(s.periodic(0), False) self.assertEqual(s.periodic(1), False) self.assertEqual(s.dimension, 3) self.assertEqual(s.rational, True) self.assertEqual(s.order(), (4, 3)) self.assertAlmostEqual(len(s.knots(0, True)), 12) self.assertAlmostEqual(len(s.knots(1, True)), 10) self.assertEqual(s.bases[1].continuity(.25), 0) self.assertEqual(s.bases[1].continuity(.75), 1) self.assertEqual(s.bases[0].continuity(.2), 2) self.assertEqual(s.bases[0].continuity(.9), 1)
def test_const_par_crv(self): # more or less random 2D surface with p=[2,2] and n=[4,3] controlpoints = [[0, 0], [-1, 1], [0, 2], [1, -1], [1, 0], [1, 1], [2, 1], [2, 2], [2, 3], [3, 0], [4, 1], [3, 2]] basis1 = BSplineBasis(3, [0, 0, 0, .4, 1, 1, 1]) basis2 = BSplineBasis(3, [0, 0, 0, 1, 1, 1]) surf = Surface(basis1, basis2, controlpoints) surf.refine(1) # try general knot in u-direction crv = surf.const_par_curve(0.3, 'u') v = np.linspace(0, 1, 13) self.assertTrue(np.allclose(surf(0.3, v), crv(v))) # try existing knot in u-direction crv = surf.const_par_curve(0.4, 'u') v = np.linspace(0, 1, 13) self.assertTrue(np.allclose(surf(0.4, v), crv(v))) # try general knot in v-direction crv = surf.const_par_curve(0.3, 'v') u = np.linspace(0, 1, 13) self.assertTrue(np.allclose(surf(u, 0.3).reshape(13, 2), crv(u))) # try start-point crv = surf.const_par_curve(0.0, 'v') u = np.linspace(0, 1, 13) self.assertTrue(np.allclose(surf(u, 0.0).reshape(13, 2), crv(u))) # try end-point crv = surf.const_par_curve(1.0, 'v') u = np.linspace(0, 1, 13) self.assertTrue(np.allclose(surf(u, 1.0).reshape(13, 2), crv(u)))
def disc(r=1, center=(0, 0, 0), normal=(0, 0, 1), type='radial', xaxis=(1, 0, 0)): """ Create a circular disc. The *type* parameter distinguishes between different parametrizations. :param float r: Radius :param array-like center: local origin :param array-like normal: local normal :param string type: The type of parametrization ('radial' or 'square') :param array-like xaxis: direction of sem, i.e. parametric start point v=0 :return: The disc :rtype: Surface """ if type == 'radial': c1 = CurveFactory.circle(r, center=center, normal=normal, xaxis=xaxis) c2 = flip_and_move_plane_geometry(c1 * 0, center, normal) result = edge_curves(c2, c1) result.swap() result.reparam((0, r), (0, 2 * pi)) return result elif type == 'square': w = 1 / sqrt(2) cp = [[-r * w, -r * w, 1], [0, -r, w], [r * w, -r * w, 1], [-r, 0, w], [0, 0, 1], [r, 0, w], [-r * w, r * w, 1], [0, r, w], [r * w, r * w, 1]] basis1 = BSplineBasis(3) basis2 = BSplineBasis(3) result = Surface(basis1, basis2, cp, True) return flip_and_move_plane_geometry(result, center, normal) else: raise ValueError('invalid type argument')
def test_force_rational(self): # more or less random 3D volume with p=[3,2,1] and n=[4,3,2] controlpoints = [[0, 0, 1], [-1, 1, 1], [0, 2, 1], [1, -1, 2], [1, 0, 2], [1, 1, 2], [2, 1, 2], [2, 2, 2], [2, 3, 2], [3, 0, 0], [4, 1, 0], [3, 2, 0], [0, 0, 3], [-1, 1, 3], [0, 2, 3], [1, -1, 5], [1, 0, 5], [1, 1, 5], [2, 1, 4], [2, 2, 4], [2, 3, 4], [3, 0, 2], [4, 1, 2], [3, 2, 2]] basis1 = BSplineBasis(4, [0, 0, 0, 0, 2, 2, 2, 2]) basis2 = BSplineBasis(3, [0, 0, 0, 1, 1, 1]) basis3 = BSplineBasis(2, [0, 0, 1, 1]) vol = Volume(basis1, basis2, basis3, controlpoints) evaluation_point1 = vol(0.23, .66, .32) control_point1 = vol[0] vol.force_rational() evaluation_point2 = vol(0.23, .66, .32) control_point2 = vol[0] # ensure that volume has not chcanged, by comparing evaluation of it self.assertAlmostEqual(evaluation_point1[0], evaluation_point2[0]) self.assertAlmostEqual(evaluation_point1[1], evaluation_point2[1]) self.assertAlmostEqual(evaluation_point1[2], evaluation_point2[2]) # ensure that we include rational weights of 1 self.assertEqual(len(control_point1), 3) self.assertEqual(len(control_point2), 4) self.assertEqual(control_point2[3], 1) self.assertEqual(vol.rational, True)
def test_evaluate(self): # creating the identity mapping by different size for all directions vol = Volume(BSplineBasis(7), BSplineBasis(6), BSplineBasis(5)) # call evaluation at a 2x3x4 grid of points u_val = np.linspace(0, 1, 2) v_val = np.linspace(0, 1, 3) w_val = np.linspace(0, 1, 4) value = vol(u_val, v_val, w_val) self.assertEqual(value.shape[0], 2) # 2 u-evaluation points self.assertEqual(value.shape[1], 3) # 3 v-evaluation points self.assertEqual(value.shape[2], 4) # 4 w-evaluation points self.assertEqual(value.shape[3], 3) # 3 dimensions (x,y,z) self.assertEqual(vol.order(), (7, 6, 5)) for i, u in enumerate(u_val): for j, v in enumerate(v_val): for k, w in enumerate(w_val): self.assertAlmostEqual(value[i, j, k, 0], u) # identity map x=u self.assertAlmostEqual(value[i, j, k, 1], v) # identity map y=v self.assertAlmostEqual(value[i, j, k, 2], w) # identity map z=w # test errors and exceptions with self.assertRaises(ValueError): val = vol(-10, .5, .5) # evalaute outside parametric domain with self.assertRaises(ValueError): val = vol(+10, .3, .3) # evalaute outside parametric domain with self.assertRaises(ValueError): val = vol(.5, -10, .123) # evalaute outside parametric domain with self.assertRaises(ValueError): val = vol(.5, +10, .123) # evalaute outside parametric domain with self.assertRaises(ValueError): val = vol(.5, .2, +10) # evalaute outside parametric domain with self.assertRaises(ValueError): val = vol(.5, .2, -10) # evalaute outside parametric domain
def test_derivative(self): # knot vector [t_1, t_2, ... t_{n+p+1}] # polynomial degree p (order-1) # n basis functions N_i(t), for i=1...n # the power basis {1,t,t^2,t^3,...} can be expressed as: # 1 = sum N_i(t) # t = sum ts_i * N_i(t) # t^2 = sum t2s_i * N_i(t) # ts_i = sum_{j=i+1}^{i+p} t_j / p # t2s_i = sum_{j=i+1}^{i+p-1} sum_{k=j+1}^{i+p} t_j*t_k / (p 2) # (p 2) = binomial coefficient # creating the mapping: # x(u,v) = u^2*v + u(1-v) # y(u,v) = v controlpoints = [[0, 0], [1.0 / 4, 0], [3.0 / 4, 0], [.75, 0], [0, 1], [0, 1], [.5, 1], [1, 1]] basis1 = BSplineBasis(3, [0, 0, 0, .5, 1, 1, 1]) basis2 = BSplineBasis(2, [0, 0, 1, 1]) surf = Surface(basis1, basis2, controlpoints) # call evaluation at a 5x4 grid of points val = surf.derivative([0, .2, .5, .6, 1], [0, .2, .4, 1], d=(1, 0)) self.assertEqual(len(val.shape), 3) # result should be wrapped in 3-index tensor self.assertEqual(val.shape[0], 5) # 5 evaluation points in u-direction self.assertEqual(val.shape[1], 4) # 4 evaluation points in v-direction self.assertEqual(val.shape[2], 2) # 2 coordinates (x,y) self.assertAlmostEqual(surf.derivative(.2, .2, d=(1, 0))[0], .88) # dx/du=2uv+(1-v) self.assertAlmostEqual(surf.derivative(.2, .2, d=(1, 0))[1], 0) # dy/du=0 self.assertAlmostEqual(surf.derivative(.2, .2, d=(0, 1))[0], -.16) # dx/dv=u^2-u self.assertAlmostEqual(surf.derivative(.2, .2, d=(0, 1))[1], 1) # dy/dv=1 self.assertAlmostEqual(surf.derivative(.2, .2, d=(1, 1))[0], -.60) # d2x/dudv=2u-1 self.assertAlmostEqual(surf.derivative(.2, .2, d=(2, 0))[0], 0.40) # d2x/dudu=2v self.assertAlmostEqual(surf.derivative(.2, .2, d=(3, 0))[0], 0.00) # d3x/du3=0 self.assertAlmostEqual(surf.derivative(.2, .2, d=(0, 2))[0], 0.00) # d2y/dv2=0 # test errors and exceptions with self.assertRaises(ValueError): val = surf.derivative(-10, .5) # evalaute outside parametric domain with self.assertRaises(ValueError): val = surf.derivative(+10, .3) # evalaute outside parametric domain with self.assertRaises(ValueError): val = surf.derivative(.5, -10) # evalaute outside parametric domain with self.assertRaises(ValueError): val = surf.derivative(.5, +10) # evalaute outside parametric domain
def test_normal(self): surf = sf.sphere(1) surf.swap() u = np.linspace(surf.start(0) + 1e-3, surf.end(0) - 1e-3, 9) v = np.linspace(surf.start(1) + 1e-3, surf.end(1) - 1e-3, 9) xpts = surf(u, v, tensor=False) npts = surf.normal(u, v, tensor=False) self.assertEqual(npts.shape, (9, 3)) # check that the normal is pointing out of the unit ball on a 9x9 evaluation grid for (x, n) in zip(xpts, npts): self.assertAlmostEqual(n[0], x[0]) self.assertAlmostEqual(n[1], x[1]) self.assertAlmostEqual(n[2], x[2]) xpts = surf(u, v) npts = surf.normal(u, v) self.assertEqual(npts.shape, (9, 9, 3)) # check that the normal is pointing out of the unit ball on a 9x9 evaluation grid for (i, j) in zip(xpts, npts): for (x, n) in zip(i, j): self.assertAlmostEqual(n[0], x[0]) self.assertAlmostEqual(n[1], x[1]) self.assertAlmostEqual(n[2], x[2]) # check single value input n = surf.normal(0, 0) self.assertEqual(len(n), 3) # test 2D surface s = Surface() n = s.normal(.5, .5) self.assertEqual(len(n), 3) self.assertAlmostEqual(n[0], 0.0) self.assertAlmostEqual(n[1], 0.0) self.assertAlmostEqual(n[2], 1.0) n = s.normal([.25, .5], [.1, .2, .3, .4, .5, .6, .7, .8, .9]) self.assertEqual(n.shape[0], 2) self.assertEqual(n.shape[1], 9) self.assertEqual(n.shape[2], 3) self.assertAlmostEqual(n[1, 4, 0], 0.0) self.assertAlmostEqual(n[1, 4, 1], 0.0) self.assertAlmostEqual(n[1, 4, 2], 1.0) n = s.normal([.25, .5, .75], [.3, .5, .9], tensor=False) self.assertEqual(n.shape, (3, 3)) for i in range(3): for j in range(2): self.assertAlmostEqual(n[i, j], 0.0) self.assertAlmostEqual(n[i, 2], 1.0) # test errors s = Surface(BSplineBasis(3), BSplineBasis(3), [[0]] * 9) # 1D-surface with self.assertRaises(RuntimeError): s.normal(.5, .5)
def get_mixed_cont_mesh(self): # Create mixed discontinuity mesh: C^0, C^0, C^{-1} nx, ny, nz = self.n Xz = self.raw.get_discontinuous_z() b1 = BSplineBasis(2, sorted(list(range(self.n[0]+1))+[0,self.n[0]])) b2 = BSplineBasis(2, sorted(list(range(self.n[1]+1))+[0,self.n[1]])) b3 = BSplineBasis(2, sorted(list(range(self.n[2]+1))*2)) true_vol = Volume(b1, b2, b3, Xz, raw=True) return true_vol
def get_cm1_mesh(self): # Create the C^{-1} mesh nx, ny, nz = self.n Xm1 = self.raw.get_discontinuous_all() b1 = BSplineBasis(2, sorted(list(range(self.n[0]+1))*2)) b2 = BSplineBasis(2, sorted(list(range(self.n[1]+1))*2)) b3 = BSplineBasis(2, sorted(list(range(self.n[2]+1))*2)) discont_vol = Volume(b1, b2, b3, Xm1) return discont_vol
def get_c0_mesh(self): # Create the C0-mesh nx, ny, nz = self.n X = self.raw.get_c0_avg() b1 = BSplineBasis(2, [0] + [i/nx for i in range(nx+1)] + [1]) b2 = BSplineBasis(2, [0] + [i/ny for i in range(ny+1)] + [1]) b3 = BSplineBasis(2, [0] + [i/nz for i in range(nz+1)] + [1]) c0_vol = volume_factory.interpolate(X, [b1, b2, b3]) return c0_vol
def loft(*curves): if len(curves) == 1: curves = curves[0] # clone input, so we don't change those references # make sure everything has the same dimension since we need to compute length curves = [c.clone().set_dimension(3) for c in curves] if len(curves)==2: return edge_curves(curves) elif len(curves)==3: # can't do cubic spline interpolation, so we'll do quadratic basis2 = BSplineBasis(3) dist = basis2.greville() else: x = [c.center() for c in curves] # create knot vector from the euclidian length between the curves dist = [0] for (x1,x0) in zip(x[1:],x[:-1]): dist.append(dist[-1] + np.linalg.norm(x1-x0)) # using "free" boundary condition by setting N'''(u) continuous at second to last and second knot knot = [dist[0]]*4 + dist[2:-2] + [dist[-1]]*4 basis2 = BSplineBasis(4, knot) n = len(curves) for i in range(n): for j in range(i+1,n): Curve.make_splines_identical(curves[i], curves[j]) basis1 = curves[0].bases[0] m = basis1.num_functions() u = basis1.greville() # parametric interpolation points v = dist # parametric interpolation points # compute matrices Nu = basis1(u) Nv = basis2(v) Nu_inv = np.linalg.inv(Nu) Nv_inv = np.linalg.inv(Nv) # compute interpolation points in physical space x = np.zeros((m,n, curves[0][0].size)) for i in range(n): x[:,i,:] = Nu * curves[i].controlpoints # solve interpolation problem cp = np.tensordot(Nv_inv, x, axes=(1,1)) cp = np.tensordot(Nu_inv, cp, axes=(1,1)) # re-order controlpoints so they match up with Surface constructor cp = cp.transpose((1, 0, 2)) cp = cp.reshape(n*m, cp.shape[2]) return Surface(basis1, basis2, cp, curves[0].rational)
def circle_segment(theta, r=1, center=(0, 0, 0), normal=(0, 0, 1), xaxis=(1, 0, 0)): """ Create a circle segment starting parallel to the rotated x-axis. :param float theta: Angle in radians :param float r: Radius :param array-like center: circle segment center :param array-like normal: normal vector to the plane that contains circle :param array-like xaxis: direction of the parametric start point t=0 :return: A quadratic rational curve :rtype: Curve :raises ValueError: If radius is not positive :raises ValueError: If theta is not in the range *[-2pi, 2pi]* """ # error test input if abs(theta) > 2 * pi: raise ValueError('theta needs to be in range [-2pi,2pi]') if r <= 0: raise ValueError('radius needs to be positive') if theta == 2 * pi: return circle(r, center, normal) # build knot vector knot_spans = int(ceil(abs(theta) / (2 * pi / 3))) knot = [0] for i in range(knot_spans + 1): knot += [i] * 2 knot += [knot_spans] # knot vector [0,0,0,1,1,2,2,..,n,n,n] knot = np.array(knot) / float( knot[-1]) * theta # set parametric space to [0,theta] n = (knot_spans - 1) * 2 + 3 # number of control points needed cp = [] t = 0 # current angle dt = float(theta) / knot_spans / 2 # angle step # build control points for i in range(n): w = 1 - (i % 2) * (1 - cos(dt) ) # weights = 1 and cos(dt) every other i x = r * cos(t) y = r * sin(t) cp += [[x, y, w]] t += dt if theta < 0: cp.reverse() result = Curve(BSplineBasis(3, np.flip(knot, 0)), cp, True) else: result = Curve(BSplineBasis(3, knot), cp, True) result.rotate(rotate_local_x_axis(xaxis, normal)) return flip_and_move_plane_geometry(result, center, normal)
def image_height(filename, N=[30, 30], p=[4, 4]): """Generate a B-spline surface approximation given by the heightmap in a grayscale image. :param str filename: Name of image file to read :param (int,int) N: Number of controlpoints in u-direction :param (int,int) p: Polynomial order (degree+1) :return: Normalized (all values between 0 and 1) heightmap approximation :rtype: :class:`splipy.Surface` """ import cv2 im = cv2.imread(filename) width = len(im) height = len(im[0]) # initialize image holder imGrey = np.zeros((len(im), len(im[0])), np.uint8) # convert to greyscale image if cv2.__version__[0] == '2': warnings.warn( FutureWarning( 'openCV v.2 will eventually be discontinued. Please update your version: \"pip install opencv-python --upgrade\"' )) cv2.cvtColor(im, cv2.cv.CV_RGB2GRAY, imGrey) else: cv2.cvtColor(im, cv2.COLOR_RGB2GRAY, imGrey) pts = [] # guess uniform evaluation points and knot vectors u = range(width) v = range(height) knot1 = [0] * 3 + range(N[0] - p[0] + 2) + [N[0] - p[0] + 1] * 3 knot2 = [0] * 3 + range(N[0] - p[0] + 2) + [N[0] - p[0] + 1] * 3 # normalize all values to be in range [0, 1] u = [float(i) / u[-1] for i in u] v = [float(i) / v[-1] for i in v] knot1 = [float(i) / knot1[-1] for i in knot1] knot2 = [float(i) / knot2[-1] for i in knot2] for j in range(height): for i in range(width): pts.append( [v[j], u[i], float(imGrey[width - i - 1][j]) / 255.0 * 1.0]) basis1 = BSplineBasis(4, knot1) basis2 = BSplineBasis(4, knot2) return surface_factory.least_square_fit(pts, [basis1, basis2], [u, v])
def test_evaluate(self): # knot vector [t_1, t_2, ... t_{n+p+1}] # polynomial degree p (order-1) # n basis functions N_i(t), for i=1...n # the power basis {1,t,t^2,t^3,...} can be expressed as: # 1 = sum N_i(t) # t = sum ts_i * N_i(t) # t^2 = sum t2s_i * N_i(t) # ts_i = sum_{j=i+1}^{i+p} t_j / p # t2s_i = sum_{j=i+1}^{i+p-1} sum_{k=j+1}^{i+p} t_j*t_k / (p 2) # (p 2) = binomial coefficient # creating the mapping: # x(u,v) = u^2*v + u(1-v) # y(u,v) = v controlpoints = [[0, 0], [1.0 / 4, 0], [3.0 / 4, 0], [.75, 0], [0, 1], [0, 1], [.5, 1], [1, 1]] basis1 = BSplineBasis(3, [0, 0, 0, .5, 1, 1, 1]) basis2 = BSplineBasis(2, [0, 0, 1, 1]) surf = Surface(basis1, basis2, controlpoints) # call evaluation at a 5x4 grid of points val = surf([0, .2, .5, .6, 1], [0, .2, .4, 1]) self.assertEqual(len(val.shape), 3) # result should be wrapped in 3-index tensor self.assertEqual(val.shape[0], 5) # 5 evaluation points in u-direction self.assertEqual(val.shape[1], 4) # 4 evaluation points in v-direction self.assertEqual(val.shape[2], 2) # 2 coordinates (x,y) # check evaluation at (0,0) self.assertAlmostEqual(val[0][0][0], 0.0) self.assertAlmostEqual(val[0][0][1], 0.0) # check evaluation at (.2,0) self.assertAlmostEqual(val[1][0][0], 0.2) self.assertAlmostEqual(val[1][0][1], 0.0) # check evaluation at (.2,.2) self.assertAlmostEqual(val[1][1][0], 0.168) self.assertAlmostEqual(val[1][1][1], 0.2) # check evaluation at (.5,.4) self.assertAlmostEqual(val[2][2][0], 0.4) self.assertAlmostEqual(val[2][2][1], 0.4) # check evaluation at (.6,1) self.assertAlmostEqual(val[3][3][0], 0.36) self.assertAlmostEqual(val[3][3][1], 1) # test errors and exceptions with self.assertRaises(ValueError): val = surf(-10, .5) # evalaute outside parametric domain with self.assertRaises(ValueError): val = surf(+10, .3) # evalaute outside parametric domain with self.assertRaises(ValueError): val = surf(.5, -10) # evalaute outside parametric domain with self.assertRaises(ValueError): val = surf(.5, +10) # evalaute outside parametric domain
def circle(r=1, center=(0,0,0), normal=(0,0,1), type='p2C0', xaxis=(1,0,0)): """ Create a circle. :param float r: Radius :param array-like center: local origin :param array-like normal: local normal :param string type: The type of parametrization ('p2C0' or 'p4C1') :param array-like xaxis: direction of sem, i.e. parametric start point t=0 :return: A periodic, quadratic rational curve :rtype: Curve :raises ValueError: If radius is not positive """ if r <= 0: raise ValueError('radius needs to be positive') if type == 'p2C0' or type == 'C0p2': w = 1.0 / sqrt(2) controlpoints = [[1, 0, 1], [w, w, w], [0, 1, 1], [-w, w, w], [-1, 0, 1], [-w, -w, w], [0, -1, 1], [w, -w, w]] knot = np.array([-1, 0, 0, 1, 1, 2, 2, 3, 3, 4, 4, 5]) / 4.0 * 2 * pi result = Curve(BSplineBasis(3, knot, 0), controlpoints, True) elif type.lower() == 'p4c1' or type.lower() == 'c1p4': w = 2*sqrt(2)/3 a = 1.0/2/sqrt(2) b = 1.0/6 * (4*sqrt(2)-1) controlpoints = [[ 1,-a, 1], [ 1, a, 1], [ b, b, w], [ a, 1, 1], [-a, 1, 1], [-b, b, w], [-1, a, 1], [-1,-a, 1], [-b,-b, w], [-a,-1, 1], [ a,-1, 1], [ b,-b, w]] knot = np.array([ -1, -1, 0, 0, 0, 1, 1, 1, 2, 2, 2, 3, 3, 3, 4, 4, 4, 5, 5]) / 4.0 * 2 * pi result = Curve(BSplineBasis(5, knot, 1), controlpoints, True) else: raise ValueError('Unkown type: %s' %(type)) result *= r result.rotate(rotate_local_x_axis(xaxis, normal)) return flip_and_move_plane_geometry(result, center, normal)
def test_periodic_boundary_evaluation(self): # Issue #66: "From left evaluations crash on periodic basis" b = BSplineBasis(3, [0,0,0,1,2,3,4,4,4]) x = b.evaluate(t=0, from_right=False) self.assertAlmostEqual(x[0,0], 0.0) x = b.evaluate(t=0, from_right=True) self.assertAlmostEqual(x[0,0], 1.0) b = BSplineBasis(3, [-1,0,0,1,2,3,4,4,5], periodic=0) x = b.evaluate(t=0, d=1, from_right=False) self.assertAlmostEqual(x[0, 0], 2.0) self.assertAlmostEqual(x[0,-1], -2.0)
def test_constructor(self): # test 3D constructor cp = [[0, 0, 0], [1, 0, 0], [0, 1, 0], [1, 1, 0]] surf = Surface(controlpoints=cp) val = surf(0.5, 0.5) self.assertEqual(val[0], 0.5) self.assertEqual(len(surf[0]), 3) # test 2D constructor cp = [[0, 0], [1, 0], [0, 1], [1, 1]] surf2 = Surface(controlpoints=cp) val = surf2(0.5, 0.5) self.assertEqual(val[0], 0.5) self.assertEqual(len(surf2[0]), 2) # test rational 2D constructor cp = [[0, 0, 1], [1, 0, 1], [0, 1, 1], [1, 1, 1]] surf3 = Surface(controlpoints=cp, rational=True) val = surf3(0.5, 0.5) self.assertEqual(val[0], 0.5) self.assertEqual(len(surf3[0]), 3) # test rational 3D constructor cp = [[0, 0, 0, 1], [1, 0, 0, 1], [0, 1, 0, 1], [1, 1, 0, 1]] surf4 = Surface(controlpoints=cp, rational=True) val = surf4(0.5, 0.5) self.assertEqual(val[0], 0.5) self.assertEqual(len(surf4[0]), 4) # test constructor with single basis b = BSplineBasis(4) surf = Surface(b, b) surf.insert_knot(.3, 'u') # change one, but not the other self.assertEqual(len(surf.knots('u')), 3) self.assertEqual(len(surf.knots('v')), 2) # TODO: Include a default constructor specifying nothing, or just polynomial degrees, or just knot vectors. # This should create identity mappings # test errors and exceptions controlpoints = [[0, 0, 1], [1, 0, 1], [0, 1, 1], [1, 1, 1]] with self.assertRaises(ValueError): basis1 = BSplineBasis(2, [1, 1, 0, 0]) basis2 = BSplineBasis(2, [0, 0, 1, 1]) surf = Surface(basis1, basis2, controlpoints) # illegal knot vector with self.assertRaises(ValueError): basis1 = BSplineBasis(2, [0, 0, .5, 1, 1]) basis2 = BSplineBasis(2, [0, 0, 1, 1]) surf = Surface(basis1, basis2, controlpoints) # too few controlpoints
def test_swap(self): # more or less random 3D volume with p=[3,2,1] and n=[4,3,2] controlpoints = [[0, 0, 1], [-1, 1, 1], [0, 2, 1], [1, -1, 2], [1, 0, 2], [1, 1, 2], [2, 1, 2], [2, 2, 2], [2, 3, 2], [3, 0, 0], [4, 1, 0], [3, 2, 0], [0, 0, 3], [-1, 1, 3], [0, 2, 3], [1, -1, 5], [1, 0, 5], [1, 1, 5], [2, 1, 4], [2, 2, 4], [2, 3, 4], [3, 0, 2], [4, 1, 2], [3, 2, 2]] basis1 = BSplineBasis(4, [0, 0, 0, 0, 2, 2, 2, 2]) basis2 = BSplineBasis(3, [0, 0, 0, 1, 1, 1]) basis3 = BSplineBasis(2, [0, 0, 1, 1]) vol = Volume(basis1, basis2, basis3, controlpoints) evaluation_point1 = vol(0.23, .56, .12) control_point1 = vol[ 1] # this is control point i=(1,0,0), when n=(4,3,2) self.assertEqual(vol.order(), (4, 3, 2)) vol.swap(0, 1) evaluation_point2 = vol(0.56, .23, .12) control_point2 = vol[ 3] # this is control point i=(0,1,0), when n=(3,4,2) self.assertEqual(vol.order(), (3, 4, 2)) # ensure that volume has not chcanged, by comparing evaluation of it self.assertAlmostEqual(evaluation_point1[0], evaluation_point2[0]) self.assertAlmostEqual(evaluation_point1[1], evaluation_point2[1]) self.assertAlmostEqual(evaluation_point1[2], evaluation_point2[2]) # check that the control points have re-ordered themselves self.assertEqual(control_point1[0], control_point2[0]) self.assertEqual(control_point1[1], control_point2[1]) self.assertEqual(control_point1[2], control_point2[2]) vol.swap(1, 2) evaluation_point3 = vol(.56, .12, .23) control_point3 = vol[ 6] # this is control point i=(0,0,1), when n=(3,2,4) self.assertEqual(vol.order(), (3, 2, 4)) # ensure that volume has not chcanged, by comparing evaluation of it self.assertAlmostEqual(evaluation_point1[0], evaluation_point3[0]) self.assertAlmostEqual(evaluation_point1[1], evaluation_point3[1]) self.assertAlmostEqual(evaluation_point1[2], evaluation_point3[2]) # check that the control points have re-ordered themselves self.assertEqual(control_point1[0], control_point3[0]) self.assertEqual(control_point1[1], control_point3[1]) self.assertEqual(control_point1[2], control_point3[2])
def test_evaluate_nontensor(self): vol = Volume(BSplineBasis(7), BSplineBasis(7), BSplineBasis(5)) u_val = [0, 0.1, 0.9, 0.3] v_val = [0.2, 0.3, 0.9, 0.4] w_val = [0.3, 0.5, 0.5, 0.0] value = vol(u_val, v_val, w_val, tensor=False) self.assertEqual(value.shape[0], 4) self.assertEqual(value.shape[1], 3) for i, (u, v, w) in enumerate(zip(u_val, v_val, w_val)): self.assertAlmostEqual(value[i, 0], u) # identity map x=u self.assertAlmostEqual(value[i, 1], v) # identity map y=v self.assertAlmostEqual(value[i, 2], w) # identity map z=w
def image_height(filename, N=[30, 30], p=[4, 4]): """Generate a B-spline surface approximation given by the heightmap in a grayscale image. :param str filename: Name of image file to read :param (int,int) N: Number of controlpoints in u-direction :param (int,int) p: Polynomial order (degree+1) :return: Normalized (all values between 0 and 1) heightmap approximation :rtype: :class:`splipy.Surface` """ import cv2 im = cv2.imread(filename) width = len(im[0]) height = len(im) # initialize image holder imGrey = np.zeros((height, width), np.uint8) # convert to greyscale image cv2.cvtColor(im, cv2.COLOR_RGB2GRAY, imGrey) # guess uniform evaluation points and knot vectors u = list(range(width)) v = list(range(height)) knot1 = [0] * (p[0] - 1) + list( range(N[0] - p[0] + 2)) + [N[0] - p[0] + 1] * (p[0] - 1) knot2 = [0] * (p[1] - 1) + list( range(N[1] - p[1] + 2)) + [N[1] - p[1] + 1] * (p[1] - 1) # normalize all values to be in range [0, 1] u = [float(i) / u[-1] for i in u] v = [float(i) / v[-1] for i in v] knot1 = [float(i) / knot1[-1] for i in knot1] knot2 = [float(i) / knot2[-1] for i in knot2] # flip and reverse image so coordinate (0,0) is at lower-left corner imGrey = imGrey.T / 255.0 imGrey = np.flip(imGrey, axis=1) x, y = np.meshgrid(u, v, indexing='ij') pts = np.stack([x, y, imGrey], axis=2) basis1 = BSplineBasis(p[0], knot1) basis2 = BSplineBasis(p[1], knot2) return surface_factory.least_square_fit(pts, [basis1, basis2], [u, v])
def rebuild(self, p, n): """ Creates an approximation to this curve by resampling it using a uniform knot vector of order *p* with *n* control points. :param int p: Polynomial discretization order :param int n: Number of control points :return: A new approximate curve :rtype: Curve """ # establish uniform open knot vector knot = [0] * p + list(range(1, n - p + 1)) + [n - p + 1] * p basis = BSplineBasis(p, knot) # set parametric range of the new basis to be the same as the old one basis.normalize() t0 = self.bases[0].start() t1 = self.bases[0].end() basis *= (t1 - t0) basis += t0 # fetch evaluation points and solve interpolation problem t = basis.greville() N = basis.evaluate(t, sparse=True) controlpoints = splinalg.spsolve(N, self.evaluate(t)) # return new resampled curve return Curve(basis, controlpoints)