예제 #1
0
def poisson_patch(bottom, right, top, left):
    from nutils import version
    if version != '3.0':
        raise ImportError('Outdated nutils version detected, only v3.0 supported. Upgrade by \"pip install --upgrade nutils\"')

    from nutils import mesh, function as fn
    from nutils import _, log

    # error test input
    if left.rational or right.rational or top.rational or bottom.rational:
        raise RuntimeError('poisson_patch not supported for rational splines')

    # these are given as a oriented loop, so make all run in positive parametric direction
    top.reverse()
    left.reverse()

    # in order to add spline surfaces, they need identical parametrization
    Curve.make_splines_identical(top, bottom)
    Curve.make_splines_identical(left, right)

    # create computational (nutils) mesh
    p1 = bottom.order(0)
    p2 = left.order(0)
    n1 = len(bottom)
    n2 = len(left)
    dim= left.dimension
    k1 = bottom.knots(0)
    k2 = left.knots(0)
    m1 = [bottom.order(0) - bottom.continuity(k) - 1 for k in k1]
    m2 = [left.order(0)   - left.continuity(k)   - 1 for k in k2]
    domain, geom = mesh.rectilinear([k1, k2])
    basis = domain.basis('spline', [p1-1, p2-1], knotmultiplicities=[m1,m2])

    # assemble system matrix
    grad      = basis.grad(geom)
    outer     = fn.outer(grad,grad)
    integrand = outer.sum(-1)
    matrix = domain.integrate(integrand, geometry=geom, ischeme='gauss'+str(max(p1,p2)+1))

    # initialize variables
    controlpoints = np.zeros((n1,n2,dim))
    rhs           = np.zeros((n1*n2))
    constraints = np.array([[np.nan]*n2]*n1)

    # treat all dimensions independently
    for d in range(dim):
        # add boundary conditions
        constraints[ 0, :] = left[  :,d]
        constraints[-1, :] = right[ :,d]
        constraints[ :, 0] = bottom[:,d]
        constraints[ :,-1] = top[   :,d]

        # solve system
        lhs = matrix.solve(rhs, constrain=np.ndarray.flatten(constraints), solver='cg', tol=state.controlpoint_absolute_tolerance)

        # wrap results into splipy datastructures
        controlpoints[:,:,d] = np.reshape(lhs, (n1,n2), order='C')

    return Surface(bottom.bases[0], left.bases[0], controlpoints, bottom.rational, raw=True)
예제 #2
0
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)
예제 #3
0
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)
예제 #4
0
def edge_curves(*curves):
    """edge_curves(curves...)

    Create the surface defined by the region between the input curves.

    In case of four input curves, these must be given in an ordered directional
    closed loop around the resulting surface.

    :param [Curve] curves: Two or four edge curves
    :return: The enclosed surface
    :rtype: Surface
    :raises ValueError: If the length of *curves* is not two or four
    """
    if len(curves) == 1: # probably gives input as a list-like single variable
        curves = curves[0]
    if len(curves) == 2:
        crv1 = curves[0].clone()
        crv2 = curves[1].clone()
        Curve.make_splines_identical(crv1, crv2)
        (n, d) = crv1.controlpoints.shape  # d = dimension + rational

        controlpoints = np.zeros((2 * n, d))
        controlpoints[:n, :] = crv1.controlpoints
        controlpoints[n:, :] = crv2.controlpoints
        linear = BSplineBasis(2)

        return Surface(crv1.bases[0], linear, controlpoints, crv1.rational)
    elif len(curves) == 4:
        # coons patch (https://en.wikipedia.org/wiki/Coons_patch)
        bottom = curves[0]
        right  = curves[1]
        top    = curves[2].clone()
        left   = curves[3].clone()  # gonna change these two, so make copies
        top.reverse()
        left.reverse()
        # create linear interpolation between opposing sides
        s1 = edge_curves(bottom, top)
        s2 = edge_curves(left, right)
        s2.swap()
        # create (linear,linear) corner parametrization
        linear = BSplineBasis(2)
        rat = s1.rational  # using control-points from top/bottom, so need to know if these are rational
        s3 = Surface(linear, linear, [bottom[0], bottom[-1], top[0], top[-1]], rat)

        # in order to add spline surfaces, they need identical parametrization
        Surface.make_splines_identical(s1, s2)
        Surface.make_splines_identical(s1, s3)
        Surface.make_splines_identical(s2, s3)

        result = s1
        result.controlpoints += s2.controlpoints
        result.controlpoints -= s3.controlpoints
        return result
    else:
        raise ValueError('Requires two or four input curves')
예제 #5
0
def loft(*curves):
    """  Generate a surface by lofting a series of curves

    The resulting surface is interpolated at all input curves and a smooth transition
    between these curves is computed as a cubic spline interpolation in the lofting
    direction. In the case that insufficient curves are provided as input (less than 4
    curves), then a quadratic or linear interpolation is performed.

    Note that the order of input curves matter as they will be interpolated in this
    particular order. Also note that the curves need to be parametrized in the same
    direction, otherwise you will encounter self-intersecting result surface as it is
    wrapping in on itself.

    :param Curve curves:  A sequence of curves to be lofted
    :return: Lofted surface
    :rtype: Surface

    Examples:

    .. code:: python

        from splipy import curve_factory, surface_factory
        crv1 = curve_factory.circle(r=1)
        crv2 = 1.5 * crv1 + [0,0,1]
        crv3 = 1.0 * crv1 + [0,0,2]
        crv4 = 1.3 * crv1 + [0,0,4]
        srf = surface_factory.loft(crv1, crv2, crv3, crv4)

        # alternatively you can provide all input curves as a list
        all_my_curves = [crv1, crv2, crv3, crv4]
        srf = surface_factory.loft(all_my_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)
예제 #6
0
def finitestrain_patch(bottom, right, top, left):
    from nutils import version
    if int(version[0]) != 4:
        raise ImportError(
            'Mismatching nutils version detected, only version 4 supported. Upgrade by \"pip install --upgrade nutils\"'
        )

    from nutils import mesh, function
    from nutils import _, log, solver

    # error test input
    if not (left.dimension == right.dimension == top.dimension ==
            bottom.dimension == 2):
        raise RuntimeError(
            'finitestrain_patch only supported for planar (2D) geometries')
    if left.rational or right.rational or top.rational or bottom.rational:
        raise RuntimeError(
            'finitestrain_patch not supported for rational splines')

    # these are given as a oriented loop, so make all run in positive parametric direction
    top.reverse()
    left.reverse()

    # in order to add spline surfaces, they need identical parametrization
    Curve.make_splines_identical(top, bottom)
    Curve.make_splines_identical(left, right)

    # create an initial mesh (correct corners) which we will morph into the right one
    p1 = bottom.order(0)
    p2 = left.order(0)
    p = max(p1, p2)
    linear = BSplineBasis(2)
    srf = Surface(linear, linear, [bottom[0], bottom[-1], top[0], top[-1]])
    srf.raise_order(p1 - 2, p2 - 2)
    for k in bottom.knots(0, True)[p1:-p1]:
        srf.insert_knot(k, 0)
    for k in left.knots(0, True)[p2:-p2]:
        srf.insert_knot(k, 1)

    # create computational mesh
    n1 = len(bottom)
    n2 = len(left)
    dim = left.dimension
    domain, geom = mesh.rectilinear(srf.knots())
    ns = function.Namespace()
    ns.basis = domain.basis('spline',
                            degree(srf),
                            knotmultiplicities=multiplicities(srf)).vector(2)
    ns.phi = domain.basis('spline',
                          degree(srf),
                          knotmultiplicities=multiplicities(srf))
    ns.eye = np.array([[1, 0], [0, 1]])
    ns.cp = controlpoints(srf)
    ns.x_i = 'cp_ni phi_n'
    ns.lmbda = 1
    ns.mu = 1

    # add total boundary conditions
    # for hard problems these will be taken in steps and multiplied by dt every
    # time (quasi-static iterations)
    constraints = np.array([[[np.nan] * n2] * n1] * dim)
    for d in range(dim):
        constraints[d, 0, :] = (left[:, d] - srf[0, :, d])
        constraints[d, -1, :] = (right[:, d] - srf[-1, :, d])
        constraints[d, :, 0] = (bottom[:, d] - srf[:, 0, d])
        constraints[d, :, -1] = (top[:, d] - srf[:, -1, d])
    # TODO: Take a close look at the logic below

    # in order to iterate, we let t0=0 be current configuration and t1=1 our target configuration
    # if solver divergeces (too large deformation), we will try with dt=0.5. If this still
    # fails we will resort to dt=0.25 until a suitable small iterations size have been found

    # dt = 1
    # t0 = 0
    # t1 = 1
    # while t0 < 1:
    # dt = t1-t0
    n = 10
    dt = 1 / n
    for i in range(n):
        # print(' ==== Quasi-static '+str(t0*100)+'-'+str(t1*100)+' % ====')
        print(' ==== Quasi-static ' + str(i / (n - 1) * 100) + ' % ====')

        # define the non-linear finite strain problem formulation
        ns.cp = np.reshape(srf[:, :, :].swapaxes(0, 1), (n1 * n2, dim),
                           order='F')
        ns.x_i = 'cp_ni phi_n'  # geometric mapping (reference geometry)
        ns.u_i = 'basis_ki ?w_k'  # displacement (unknown coefficients w_k)
        ns.X_i = 'x_i + u_i'  # displaced geometry
        ns.strain_ij = '.5 (u_i,j + u_j,i + u_k,i u_k,j)'
        ns.stress_ij = 'lmbda strain_kk eye_ij + 2 mu strain_ij'

        # try:
        residual = domain.integral(ns.eval_n('stress_ij basis_ni,j d:X'),
                                   degree=2 * p)
        cons = np.ndarray.flatten(constraints * dt, order='C')
        lhs = solver.newton('w', residual, constrain=cons).solve(
            tol=state.controlpoint_absolute_tolerance, maxiter=8)

        # store the results on a splipy object and continue
        geom = lhs.reshape((n2, n1, dim), order='F')
        srf[:, :, :] += geom.swapaxes(0, 1)

        # t0 += dt
        # t1 = 1
        # except solver.SolverError: # newton method fail to converge, try a smaller step length 'dt'
        # t1 = (t1+t0)/2
    return srf
예제 #7
0
def elasticity_patch(bottom, right, top, left):
    from nutils import version
    if int(version[0]) != 4:
        raise ImportError(
            'Mismatching nutils version detected, only version 4 supported. Upgrade by \"pip install --upgrade nutils\"'
        )

    from nutils import mesh, function, solver
    from nutils import _, log

    # error test input
    if not (left.dimension == right.dimension == top.dimension ==
            bottom.dimension == 2):
        raise RuntimeError(
            'elasticity_patch only supported for planar (2D) geometries')
    if left.rational or right.rational or top.rational or bottom.rational:
        raise RuntimeError(
            'elasticity_patch not supported for rational splines')

    # these are given as a oriented loop, so make all run in positive parametric direction
    top.reverse()
    left.reverse()

    # in order to add spline surfaces, they need identical parametrization
    Curve.make_splines_identical(top, bottom)
    Curve.make_splines_identical(left, right)

    # create computational (nutils) mesh
    p1 = bottom.order(0)
    p2 = left.order(0)
    n1 = len(bottom)
    n2 = len(left)
    dim = left.dimension
    k1 = bottom.knots(0)
    k2 = left.knots(0)
    m1 = [bottom.order(0) - bottom.continuity(k) - 1 for k in k1]
    m2 = [left.order(0) - left.continuity(k) - 1 for k in k2]
    domain, geom = mesh.rectilinear([k1, k2])

    # assemble system matrix
    ns = function.Namespace()
    ns.x = geom
    ns.basis = domain.basis('spline',
                            degree=[p1 - 1, p2 - 1],
                            knotmultiplicities=[m1, m2]).vector(2)
    ns.eye = np.array([[1, 0], [0, 1]])
    ns.lmbda = 1
    ns.mu = .3
    # ns.u_i = 'basis_ni ?lhs_n'
    # ns.strain_ij = '(u_i,j + u_j,i) / 2'
    # ns.stress_ij = 'lmbda strain_kk eye_ij + 2 mu strain_ij'
    ns.strain_nij = '(basis_ni,j + basis_nj,i) / 2'
    ns.stress_nij = 'lmbda strain_nkk eye_ij + 2 mu strain_nij'

    # construct matrix and right hand-side
    matrix = domain.integrate(ns.eval_nm('strain_nij stress_mij d:x'),
                              ischeme='gauss' + str(max(p1, p2) + 1))
    rhs = np.zeros((n1 * n2 * dim))

    # add boundary conditions
    constraints = np.array([[[np.nan] * n2] * n1] * dim)
    for d in range(dim):
        constraints[d, 0, :] = left[:, d]
        constraints[d, -1, :] = right[:, d]
        constraints[d, :, 0] = bottom[:, d]
        constraints[d, :, -1] = top[:, d]

    # solve system
    lhs = matrix.solve(rhs,
                       constrain=np.ndarray.flatten(constraints, order='C'))

    # rewrap results into splipy datastructures
    controlpoints = np.reshape(lhs, (dim, n1, n2), order='C')
    controlpoints = controlpoints.swapaxes(0, 2)
    controlpoints = controlpoints.swapaxes(0, 1)

    return Surface(bottom.bases[0],
                   left.bases[0],
                   controlpoints,
                   bottom.rational,
                   raw=True)
예제 #8
0
def edge_curves(*curves, **kwargs):
    """  Create the surface defined by the region between the input curves.

    In case of four input curves, these must be given in an ordered directional
    closed loop around the resulting surface.

    :param [Curve] curves: Two or four edge curves
    :param string type: The method used for interior computation ('coons', 'poisson', 'elasticity' or 'finitestrain')
    :return: The enclosed surface
    :rtype: Surface
    :raises ValueError: If the length of *curves* is not two or four
    """
    type = kwargs.get('type', 'coons')
    if len(curves) == 1:  # probably gives input as a list-like single variable
        curves = curves[0]
    if len(curves) == 2:
        crv1 = curves[0].clone()
        crv2 = curves[1].clone()
        Curve.make_splines_identical(crv1, crv2)
        (n, d) = crv1.controlpoints.shape  # d = dimension + rational

        controlpoints = np.zeros((2 * n, d))
        controlpoints[:n, :] = crv1.controlpoints
        controlpoints[n:, :] = crv2.controlpoints
        linear = BSplineBasis(2)

        return Surface(crv1.bases[0], linear, controlpoints, crv1.rational)
    elif len(curves) == 4:
        # reorganize input curves so they form a directed loop around surface
        rtol = state.controlpoint_relative_tolerance
        atol = state.controlpoint_absolute_tolerance
        mycurves = [c.clone() for c in curves
                    ]  # wrap into list and clone all since we're changing them
        dim = np.max([c.dimension for c in mycurves])
        rat = np.any([c.rational for c in mycurves])
        for i in range(4):
            for j in range(i + 1, 4):
                Curve.make_splines_compatible(mycurves[i], mycurves[j])

        if not (np.allclose(
                mycurves[0][-1], mycurves[1][0], rtol=rtol, atol=atol)
                and np.allclose(
                    mycurves[1][-1], mycurves[2][0], rtol=rtol, atol=atol)
                and np.allclose(
                    mycurves[2][-1], mycurves[3][0], rtol=rtol, atol=atol)
                and np.allclose(
                    mycurves[3][-1], mycurves[0][0], rtol=rtol, atol=atol)):
            reorder = [mycurves[0]]
            del mycurves[0]
            for j in range(3):
                found_match = False
                for i in range(len(mycurves)):
                    if (np.allclose(reorder[j][-1],
                                    mycurves[i][0],
                                    rtol=rtol,
                                    atol=atol)):
                        reorder.append(mycurves[i])
                        del mycurves[i]
                        found_match = True
                        break
                    elif (np.allclose(reorder[j][-1],
                                      mycurves[i][-1],
                                      rtol=rtol,
                                      atol=atol)):
                        reorder.append(mycurves[i].reverse())
                        del mycurves[i]
                        found_match = True
                        break
                if not found_match:
                    raise RuntimeError(
                        'Curves do not form a closed loop (end-points do not match)'
                    )
            mycurves = reorder
        if type == 'coons':
            return coons_patch(*mycurves)
        elif type == 'poisson':
            return poisson_patch(*mycurves)
        elif type == 'elasticity':
            return elasticity_patch(*mycurves)
        elif type == 'finitestrain':
            return finitestrain_patch(*mycurves)
        else:
            raise ValueError('Unknown type parameter')
    else:
        raise ValueError('Requires two or four input curves')