예제 #1
0
def load_mesh(filename, extension_hint=None, drop_zero_dim=False):
    """ Load mesh from a file.

    Args:
        filename: Input filename.  File format is auto detected based on
            extension.
        drop_zero_dim (bool): If true, convert flat 3D mesh into 2D mesh.

    Returns:
        A :py:class:`Mesh` object representing the loaded mesh.
    """

    ext = os.path.splitext(filename)[1] if extension_hint is None else extension_hint

    if ext == ".geogram":
        return Mesh(PyMesh.load_geogram_mesh(filename));
    if not os.path.exists(filename):
        raise IOError("File not found: {}".format(filename));
    factory = PyMesh.MeshFactory();
    if extension_hint is None:
        factory.load_file(filename);
    else:
        factory.load_file_with_hint(filename, extension_hint)
    if drop_zero_dim:
        factory.drop_zero_dim();
    return Mesh(factory.create());
예제 #2
0
 def __init__(self, cell_size, dim=3):
     self.dim = dim
     if dim == 3:
         self.raw_grid = PyMesh.VoxelGrid3D(cell_size)
     elif dim == 2:
         self.raw_grid = PyMesh.VoxelGrid2D(cell_size)
     else:
         raise NotImplementedError("Unsupported dim: {}".format(dim))
예제 #3
0
def form_mesh(vertices, faces, voxels=None):
    """ Convert raw mesh data into a Mesh object.

    Args:
        vertices (`numpy.ndarray`): ndarray of floats with size (num_vertices,
            dim).
        faces: ndarray of ints with size (num_faces, vertex_per_face).
        voxels: optional ndarray of ints with size (num_voxels,
            vertex_per_voxel).  Use ``None`` for forming surface meshes.

    Returns:
        A :py:class:`Mesh` object formed by the inputs.

    Example:

        >>> vertices = np.array([
        ...     [0.0, 0.0],
        ...     [1.0, 0.0],
        ...     [1.0, 1.0],
        ...     [0.0, 1.0],
        ...     ]);
        >>> faces = np.array([
        ...     [0, 1, 2],
        ...     [0, 2, 3],
        ...     ]);
        >>> mesh = pymesh.form_mesh(vertices, faces);
    """
    voxels = deduce_voxel_type(faces, voxels)
    faces = deduce_face_type(faces, voxels)

    factory = PyMesh.MeshFactory()
    factory.load_data(vertices.ravel(order="C"), faces.ravel(order="C"),
                      voxels.ravel(order="C"), vertices.shape[1],
                      faces.shape[1], voxels.shape[1])
    return Mesh(factory.create())
예제 #4
0
def snap_rounding(points, segments, pixel_size, use_iterative=True):
    """ 2D snap rounding.

    Args:
        points (``numpy.ndarray``): Input points.
        segments (``numpy.ndarray``): Input segments.
        pixel_size (``float``): Pixel size.
        use_iterative (``bool``): Whether to use iterative snap rounding.

    Returns:
        2 values are returned.

            * ``vertices``: Snap rounded vertices.
            * ``edges``: Edge connecting vertices.
    """
    engine = PyMesh.SnapRounding2()
    engine.points = points
    engine.segments = segments
    engine.run(pixel_size, use_iterative)

    vertices = engine.vertices
    edges = engine.edges

    vertices, edges, __ = remove_duplicated_vertices_raw(vertices,
                                                         edges,
                                                         tol=pixel_size / 2)
    return vertices, edges
예제 #5
0
 def _extra_info(self):
     try:
         return self.__extra_info
     except AttributeError:
         self.__extra_info = PyMesh.MeshChecker(
                 self.vertices, self.faces, self.voxels)
         return self.__extra_info
예제 #6
0
def form_mesh(vertices, faces, voxels=None):
    """ Convert raw mesh data into a Mesh object.

    Args:
        vertices (`numpy.ndarray`): ndarray of floats with size (num_vertices,
            dim).
        faces: ndarray of ints with size (num_faces, vertex_per_face).
        voxels: optional ndarray of ints with size (num_voxels,
            vertex_per_voxel).  Use ``None`` for forming surface meshes.

    Returns:
        A Mesh object.
    """
    voxels = deduce_voxel_type(faces, voxels);
    faces = deduce_face_type(faces, voxels);

    factory = PyMesh.MeshFactory();
    factory.load_data(
            vertices.ravel(order="C"),
            faces.ravel(order="C"),
            voxels.ravel(order="C"),
            vertices.shape[1],
            faces.shape[1],
            voxels.shape[1]);
    return Mesh(factory.create());
예제 #7
0
파일: meshio.py 프로젝트: qnzhou/PyMesh
def load_mesh(filename, extension_hint=None, drop_zero_dim=False):
    """ Load mesh from a file.

    Args:
        filename: Input filename.  File format is auto detected based on
            extension.
        drop_zero_dim (bool): If true, convert flat 3D mesh into 2D mesh.

    Returns:
        A :py:class:`Mesh` object representing the loaded mesh.
    """

    ext = os.path.splitext(filename)[1] if extension_hint is None else extension_hint

    if ext == ".geogram":
        return Mesh(PyMesh.load_geogram_mesh(filename));
    if not os.path.exists(filename):
        raise IOError("File not found: {}".format(filename));
    factory = PyMesh.MeshFactory();
    if extension_hint is None:
        factory.load_file(filename);
    else:
        factory.load_file_with_hint(filename, extension_hint)
    if drop_zero_dim:
        factory.drop_zero_dim();
    return Mesh(factory.create());
예제 #8
0
def save_mesh(filename, mesh, *attributes, **setting):
    """ Save mesh to file.

    Args:
        filename (:py:class:`str`): Output file.  File format is auto detected from extension.
        mesh (:py:class:`Mesh`): Mesh object.
        *attributes (:py:class:`list`): (optional) Attribute names to be saved.
            This field would be ignored if the output format does not support
            attributes (e.g. **.obj** and **.stl** files)
        **setting (:py:class:`dict`): (optional) The following keys are recognized.

            * ascii: whether to use ascii encoding, default is false.
            * use_float: store scalars as float instead of double, default is
              false.
            * anonymous: whether to indicate the file is generated by PyMesh.

    Raises:
        KeyError: Attributes cannot be found in mesh.

    Example:

        >>> box = pymesh.generate_box_mesh();
        >>> pymesh.save_mesh("tmp.stl", box, ascii=True);
    """
    ext = os.path.splitext(filename)[1];
    if ext == ".geogram":
        PyMesh.save_geogram_mesh(filename, mesh.raw_mesh);
        return;
    elif ext == ".svg":
        save_svg(filename, mesh);
        return;
    writer = PyMesh.MeshWriter.create(filename);
    for attr in attributes:
        if not mesh.has_attribute(attr):
            raise KeyError("Attribute {} is not found in mesh".format(attr));
        writer.with_attribute(attr);
    if setting.get("ascii", False):
        writer.in_ascii();
    if setting.get("use_float", False):
        writer.use_float();
    if setting.get("anonymous", False):
        writer.set_anonymous();
    writer.use_default_physical_tag( 
        setting.get("use_default_physical_tag", False )
    )
    writer.write_mesh(mesh.raw_mesh);
예제 #9
0
def is_vertex_manifold(mesh):
    """ Vertex manifold check.  The result is stored as a per-vertex scalar
    field named "vertex_manifold".
    """

    vertex_manifold = PyMesh.is_vertex_manifold(mesh.faces);
    mesh.add_attribute("vertex_manifold");
    mesh.set_attribute("vertex_manifold", vertex_manifold);
예제 #10
0
def is_edge_manifold(mesh):
    """ Edge manifold check.  The result is stored as a per-edge scalar
    field named "edge_manifold".
    """

    edge_manifold = PyMesh.is_edge_manifold(mesh.faces)
    mesh.add_attribute("edge_manifold")
    mesh.set_attribute("edge_manifold", edge_manifold)
예제 #11
0
def is_delaunay(mesh):
    """ A thin wrapper of ``is_delaunay_raw``.
    """
    if mesh.num_voxels == 0:
        return np.zeros(0)
    if mesh.vertex_per_voxel != 4:
        raise NotImplementedError("Delaunay property computation expect a tet mesh.")
    return PyMesh.is_delaunay(mesh.vertices, mesh.voxels)
예제 #12
0
def get_tet_orientations(mesh):
    """ A thin wrapper of ``get_tet_orientations_raw``.
    """
    if mesh.num_voxels == 0:
        return np.zeros(0);
    if mesh.vertex_per_voxel != 4:
        raise NotImplementedError("Distortion computation expect a tet mesh.");
    return PyMesh.get_tet_orientations(mesh.vertices, mesh.voxels);
예제 #13
0
 def tile_with_guide_mesh(self, mesh, parameters=None):
     self.__check_base_pattern()
     tiler = PyMesh.WireTiler(self.raw_pattern)
     if parameters is None:
         parameters = Parameters(self.pattern)
     tiler.with_parameters(parameters.raw_parameters)
     self.raw_wire_network = tiler.tile_with_guide_mesh(mesh.raw_mesh)
     self.__apply_vertex_offset()
예제 #14
0
def is_vertex_manifold(mesh):
    """ Vertex manifold check.  The result is stored as a per-vertex scalar
    field named "vertex_manifold".
    """

    vertex_manifold = PyMesh.is_vertex_manifold(mesh.faces)
    mesh.add_attribute("vertex_manifold")
    mesh.set_attribute("vertex_manifold", vertex_manifold)
예제 #15
0
def get_tet_orientations(mesh):
    """ A thin wrapper of ``get_tet_orientations_raw``.
    """
    if mesh.num_voxels == 0:
        return np.zeros(0)
    if mesh.vertex_per_voxel != 4:
        raise NotImplementedError("Orientation computation expect a tet mesh.")
    return PyMesh.get_tet_orientations(mesh.vertices, mesh.voxels)
예제 #16
0
    def load_mesh(self, mesh_file):
        mesh_file = os.path.join(self.data_dir, mesh_file)
        if not os.path.exists(mesh_file):
            raise IOError("mesh file {} does not exist!".format(mesh_file))

        factory = PyMesh.MeshFactory()
        factory.load_file(mesh_file)
        return factory.create_shared()
예제 #17
0
def is_edge_manifold(mesh):
    """ Edge manifold check.  The result is stored as a per-edge scalar
    field named "edge_manifold".
    """

    edge_manifold = PyMesh.is_edge_manifold(mesh.faces);
    mesh.add_attribute("edge_manifold");
    mesh.set_attribute("edge_manifold", edge_manifold);
예제 #18
0
 def tile_with_guide_bbox(self, bbox_min, bbox_max, reps, parameters=None):
     self.__check_base_pattern()
     tiler = PyMesh.WireTiler(self.raw_pattern)
     if parameters is None:
         parameters = Parameters(self.pattern)
     tiler.with_parameters(parameters.raw_parameters)
     self.raw_wire_network = tiler.tile_with_guide_bbox(
         np.array(bbox_min, dtype=float), np.array(bbox_max, dtype=float),
         np.array(reps, dtype=int))
     self.__apply_vertex_offset()
예제 #19
0
파일: meshio.py 프로젝트: qnzhou/PyMesh
def save_mesh(filename, mesh, *attributes, **setting):
    """ Save mesh to file.

    Args:
        filename (:py:class:`str`): Output file.  File format is auto detected from extension.
        mesh (:py:class:`Mesh`): Mesh object.
        *attributes (:py:class:`list`): (optional) Attribute names to be saved.
            This field would be ignored if the output format does not support
            attributes (e.g. **.obj** and **.stl** files)
        **setting (:py:class:`dict`): (optional) The following keys are recognized.

            * ascii: whether to use ascii encoding, default is false.
            * use_float: store scalars as float instead of double, default is
              false.
            * anonymous: whether to indicate the file is generated by PyMesh.

    Raises:
        KeyError: Attributes cannot be found in mesh.

    Example:

        >>> box = pymesh.generate_box_mesh();
        >>> pymesh.save_mesh("tmp.stl", box, ascii=True);
    """
    ext = os.path.splitext(filename)[1];
    if ext == ".geogram":
        PyMesh.save_geogram_mesh(filename, mesh.raw_mesh);
        return;
    elif ext == ".svg":
        save_svg(filename, mesh);
        return;
    writer = PyMesh.MeshWriter.create(filename);
    for attr in attributes:
        if not mesh.has_attribute(attr):
            raise KeyError("Attribute {} is not found in mesh".format(attr));
        writer.with_attribute(attr);
    if setting.get("ascii", False):
        writer.in_ascii();
    if setting.get("use_float", False):
        writer.use_float();
    if setting.get("anonymous", False):
        writer.set_anonymous();
    writer.write_mesh(mesh.raw_mesh);
예제 #20
0
 def tile_with_mixed_patterns(self,
                              mesh,
                              per_vertex_thickness=False,
                              isotropic_dofs=True):
     self.__check_base_patterns()
     tiler = PyMesh.WireTiler(self.raw_patterns[0])
     self.raw_wire_network = tiler.tile_with_mixed_patterns(
         self.raw_patterns, mesh.raw_mesh, per_vertex_thickness,
         isotropic_dofs)
     self.__apply_vertex_offset()
예제 #21
0
파일: face_utils.py 프로젝트: qnzhou/PyMesh
def is_colinear(v0, v1, v2):
    """ Return true if ``v0``, ``v1`` and ``v2`` are colinear.
    Colinear check is done using exact predicates.

    Args:
        v0 (``numpy.ndarray``): vector of size 2 or 3.
        v1 (``numpy.ndarray``): vector of size 2 or 3.
        v2 (``numpy.ndarray``): vector of size 2 or 3.

    Return:
        A boolean indicating whether ``v0``, ``v1`` and ``v2`` are colinear.
    """
    dim = len(v0);
    if dim == 2:
        return PyMesh.is_colinear_2D(v0, v1, v2);
    elif dim == 3:
        return PyMesh.is_colinear_3D(v0, v1, v2);
    else:
        raise NotImplementedError("Supported dimention {}".format(dim));
예제 #22
0
def is_colinear(v0, v1, v2):
    """ Return true if ``v0``, ``v1`` and ``v2`` are colinear.
    Colinear check is done using exact predicates.

    Args:
        v0 (``numpy.ndarray``): vector of size 2 or 3.
        v1 (``numpy.ndarray``): vector of size 2 or 3.
        v2 (``numpy.ndarray``): vector of size 2 or 3.

    Return:
        A boolean indicating whether ``v0``, ``v1`` and ``v2`` are colinear.
    """
    dim = len(v0)
    if dim == 2:
        return PyMesh.is_colinear_2D(v0, v1, v2)
    elif dim == 3:
        return PyMesh.is_colinear_3D(v0, v1, v2)
    else:
        raise NotImplementedError("Supported dimention {}".format(dim))
예제 #23
0
def orient_2D(p1, p2, p3):
    """ Determine the orientation 2D points p1, p2, p3

    Args:
        p1,p2,p3: 2D points.

    Returns:
        positive if (p1, p2, p3) is in counterclockwise order.
        negative if (p1, p2, p3) is in clockwise order.
        0.0 if they are collinear.
    """
    return PyMesh.orient2d(p1, p2, p3);
예제 #24
0
def get_degenerated_faces_raw(vertices, faces):
    """ Return indices of degenerated faces.
    A face is degenerated if all its 3 corners are colinear.

    Args:
        vertices (``numpy.ndarray``): Vertex matrix.
        faces (``numpy.ndarray``): Face matrix.

    Returns:
        A ``numpy.ndarray`` of indices of degenerated faces.
    """
    return np.array(PyMesh.get_degenerated_faces(vertices, faces))
예제 #25
0
def orient_3D(p1, p2, p3, p4):
    """ Determine the orientation 3D points p1, p2, p3, p4.

    Args:
        p1,p2,p3,p4: 3D points.

    Returns:
        positive if p4 is below the plane formed by (p1, p2, p3).
        negative if p4 is above the plane formed by (p1, p2, p3).
        0.0 if they are coplanar.
    """
    return PyMesh.orient3d(p1, p2, p3, p4);
예제 #26
0
파일: face_utils.py 프로젝트: qnzhou/PyMesh
def get_degenerated_faces_raw(vertices, faces):
    """ Return indices of degenerated faces.
    A face is degenerated if all its 3 corners are colinear.

    Args:
        vertices (``numpy.ndarray``): Vertex matrix.
        faces (``numpy.ndarray``): Face matrix.

    Returns:
        A ``numpy.ndarray`` of indices of degenerated faces.
    """
    return np.array(PyMesh.get_degenerated_faces(vertices, faces));
def save_matrix(filename, matrix, in_ascii=False):
    """ Save matrix into file in `.dmat`_ format.

    Args:
        filename (``str``): Output file name.
        matrix (:py:class:`numpy.ndarray`): The matrix to save.
        in_ascii (`boolean`): Whether to save matrix in ASCII.  Default is
            false, which saves in binary format to save space.

    .. _.dmat: http://libigl.github.io/libigl/file-formats/dmat.html
    """
    return PyMesh.save_matrix(filename, matrix, in_ascii)
예제 #28
0
def in_sphere(p1, p2, p3, p4, p5):
    """ Determine if p5 is in the sphere formed by p1, p2, p3, p4.

    Args:
        p1,p2,p3,p4,p5: 3D points.  ``orient_3D(p1, p2, p3, p4)`` must be
            positive, otherwise the result will be flipped.

    Returns:
        positive p5 is inside of the sphere.
        negative p5 is outside of the sphere.
        0.0 if they are cospherical.
    """
    return PyMesh.insphere(p1, p2, p3, p4, p5);
예제 #29
0
def convert_to_voxel_attribute(mesh, attr):
    """ Convert attribute ``attr`` from either per-vertex or per-face attribute
    into per-voxel attribute.

    Args:
        mesh (:py:class:`Mesh`): Input mesh.
        attr (``numpy.ndarray``): #voxel by k matrix of floats.

    Returns:
        Per-voxel attribute. The value at a voxel will be the average of
        the values at its neighboring vertices or faces.
    """
    return PyMesh.convert_to_voxel_attribute(mesh.raw_mesh, attr)
예제 #30
0
    def load_mesh(self, mesh_file):
        mesh_file = os.path.join(self.data_dir, mesh_file)
        if not os.path.exists(mesh_file):
            raise IOError("mesh file {} does not exist!".format(mesh_file))

        factory = PyMesh.MeshFactory()
        factory.load_file(mesh_file)
        factory.with_connectivity("all")
        factory.with_attribute("face_normal")
        factory.with_attribute("vertex_normal")
        factory.with_attribute("face_area")
        factory.with_attribute("voxel_volume")
        return factory.create_shared()
예제 #31
0
def convert_to_voxel_attribute(mesh, attr):
    """ Convert attribute ``attr`` from either per-vertex or per-face attribute
    into per-voxel attribute.

    Args:
        mesh (:py:class:`Mesh`): Input mesh.
        attr (``numpy.ndarray``): #voxel by k matrix of floats.

    Returns:
        Per-voxel attribute. The value at a voxel will be the average of
        the values at its neighboring vertices or faces.
    """
    return PyMesh.convert_to_voxel_attribute(mesh.raw_mesh, attr);
예제 #32
0
def in_circle(p1, p2, p3, p4):
    """ Determine if p4 is in the circle formed by p1, p2, p3.

    Args:
        p1,p2,p3,p4: 2D points.  ``orient_2D(p1, p2, p3)`` must be postive,
            otherwise the result will be flipped.

    Returns:
        positive p4 is inside of the circle.
        negative p4 is outside of the circle.
        0.0 if they are cocircular.
    """
    return PyMesh.incircle(p1, p2, p3, p4);
예제 #33
0
def get_tet_orientations_raw(vertices, tets):
    """ Compute orientation of each tet.

    Args:
        vertice (``numpy.ndarray``): n by 3 matrix representing vertices.
        tets (``numpy.ndarray``): m by 4 matrix of vertex indices representing tets.

    Returns:
        A list of m floats, where
            * Positive number => tet is positively oriented.
            *               0 => tet is degenerate.
            * Negative number => tet is inverted.
    """
    return PyMesh.get_tet_orientations(vertices, tets)
예제 #34
0
def get_tet_orientations_raw(vertices, tets):
    """ Compute orientation of each tet.

    Args:
        vertice (``numpy.ndarray``): n by 3 matrix representing vertices.
        tets (``numpy.ndarray``): m by 4 matrix of vertex indices representing tets.

    Returns:
        A list of m floats, where
            * Positive number => tet is positively oriented.
            *               0 => tet is degenerate.
            * Negative number => tet is inverted.
    """
    return PyMesh.get_tet_orientations(vertices, tets);
예제 #35
0
def is_delaunay_raw(vertices, tets):
    """ Compute whether each tet is strictly locally Delaunay, cospherical or
    not locally Delaunay.

    Args:
        vertices (``numpy.ndarray``): n by 3 matrix representing vertices.
        tets (``numpy.ndarray``): m by 4 matrix of vertex indices representing tets.

    Returns:
        A list of m floats representing result for each tet:
            * 1 => tet is strictly locally Delaunay.
            * 0 => tet is cospherical with an adjacent tet.
            * 1 => tet is not locally Delaunay.
    """
    return PyMesh.is_delaunay(vertices, tets)
예제 #36
0
def detect_self_intersection(mesh):
    """ Detect all self-intersections.

    Args:
        mesh (:class:`Mesh`): The input mesh.

    Returns:
        :py:class:`numpy.ndarray`: A n by 2 array of face indices.  Each
        row contains the indices of two intersecting faces. :math:`n` is
        the number of intersecting face pairs.
    """

    detector = PyMesh.SelfIntersection(mesh.vertices, mesh.faces)
    detector.detect_self_intersection()
    intersecting_faces = detector.get_self_intersecting_pairs()
    return intersecting_faces
예제 #37
0
파일: face_utils.py 프로젝트: qnzhou/PyMesh
def get_triangle_orientations_raw(vertices, faces):
    """ Compute orientation of each triangle.  A triangle is considered as
    positively oriented iff its three corners are counterclockwise ordered.

    Args:
        vertice (``numpy.ndarray``): n by 3 matrix representing vertices.
        faces (``numpy.ndarray``): m by 3 matrix of vertex indices representing
            triangles.

    Returns:
        A list of m floats, where
            * Positive number => triangle is positively oriented.
            *               0 => triangle is degenerate.
            * Negative number => triangle is negatively oriented.
    """
    return PyMesh.get_triangle_orientations(vertices, faces);
예제 #38
0
def get_triangle_orientations_raw(vertices, faces):
    """ Compute orientation of each triangle.  A triangle is considered as
    positively oriented iff its three corners are counterclockwise ordered.

    Args:
        vertice (``numpy.ndarray``): n by 3 matrix representing vertices.
        faces (``numpy.ndarray``): m by 3 matrix of vertex indices representing
            triangles.

    Returns:
        A list of m floats, where
            * Positive number => triangle is positively oriented.
            *               0 => triangle is degenerate.
            * Negative number => triangle is negatively oriented.
    """
    return PyMesh.get_triangle_orientations(vertices, faces)
예제 #39
0
def straight_skeleton(vertices, edges):
    """ Compute 2D straight skeleton.

    Args:
        vertices (``numpy.ndarray``): Vertices of the input boundary.
        edges (``numpy.ndarray``): Edges of the input boundary.

    Returns:
        2 values are returned:

            * ``points``: Straight skeleton vertices.
            * ``segments``: Straight skeleton segments.
    """
    engine = PyMesh.StraightSkeleton()
    engine.run(vertices, edges)
    return engine.points, engine.segments
예제 #40
0
def load_mesh(filename, drop_zero_dim=False):
    """ Load mesh from a file.

    Args:
        filename: Input filename.  File format is auto detected based on
            extension.
        drop_zero_dim (bool): If true, convert flat 3D mesh into 2D mesh.

    Returns:
        A :py:class:`Mesh` object representing the loaded mesh.
    """
    if not os.path.exists(filename):
        raise IOError("File not found: {}".format(filename))
    factory = PyMesh.MeshFactory()
    factory.load_file(filename)
    if drop_zero_dim:
        factory.drop_zero_dim()
    return Mesh(factory.create())
예제 #41
0
def convert_to_voxel_attribute_from_name(mesh, name):
    """ Same as :py:func:`convert_to_voxel_attribute` except looking up
    attribute values from the input ``mesh`` using ``name``.
    """
    return PyMesh.convert_to_voxel_attribute_from_name(mesh.raw_mesh, name);
예제 #42
0
파일: predicates.py 프로젝트: qnzhou/PyMesh
import PyMesh

"""
This module wraps the exact predicates Jonathan Richard Shewchuk.
"""

# The init function would be called when predicates module is imported.
PyMesh.exactinit();

def orient_2D(p1, p2, p3):
    """ Determine the orientation 2D points p1, p2, p3

    Args:
        p1,p2,p3: 2D points.

    Returns:
        positive if (p1, p2, p3) is in counterclockwise order.
        negative if (p1, p2, p3) is in clockwise order.
        0.0 if they are collinear.
    """
    return PyMesh.orient2d(p1, p2, p3);

def orient_3D(p1, p2, p3, p4):
    """ Determine the orientation 3D points p1, p2, p3, p4.

    Args:
        p1,p2,p3,p4: 3D points.

    Returns:
        positive if p4 is below the plane formed by (p1, p2, p3).
        negative if p4 is above the plane formed by (p1, p2, p3).
예제 #43
0
def cut_to_manifold(mesh):
    """ Cut an input mesh along nonmanifold edges to it becomes manifold.
    Note that cutting will produce duplicated vertices.
    """
    return Mesh(PyMesh.cut_to_manifold(mesh.raw_mesh));
예제 #44
0
파일: aabb_tree.py 프로젝트: qnzhou/PyMesh
def signed_distance_to_mesh(mesh, pts, engine="igl", bvh=None):
    """ Compute the signed distance from a set of points to a mesh.

    Args:
        mesh (:class:`Mesh`): A input mesh.
        pts (:class:`numpy.ndarray`): A :math:`N \\times dim` array of query
            points.
        engine (``string``): BVH engine name. Valid choices are "cgal",
            "geogram", "igl" if all dependencies are used. The default is
            "auto" where an available engine is automatically picked.
        bvh (:class:`BVH`): BVH engine instance (optional)

    Returns:
        Three values are returned.

            * ``signed_distances``: signed (unsquared) distances from each
                                    point to mesh.
            * ``face_indices``  : the closest face to each point.
            * ``closest_points``: the point on mesh that is closest to each
                                  query point.
    """

    if not bvh:
        bvh = BVH(engine, mesh.dim);
        bvh.load_mesh(mesh);

    # get face normals
    try:
        face_normals = np.reshape(mesh.get_attribute("face_normals"),
                                  np.int32(mesh.get_attribute("face_normals_shape")))
    except RuntimeError:
        face_normals = PyMesh.face_normals(mesh.vertices, mesh.faces)
        mesh.add_attribute("face_normals")
        mesh.add_attribute("face_normals_shape")
        mesh.set_attribute("face_normals", face_normals)
        mesh.set_attribute("face_normals_shape", np.array(face_normals.shape))

    # get vertex normals
    try:
        vertex_normals = np.reshape(mesh.get_attribute("vertex_normals"),
                                    np.int32(mesh.get_attribute("vertex_normals_shape")))
    except RuntimeError:
        vertex_normals = PyMesh.vertex_normals(mesh.vertices, mesh.faces, face_normals)
        mesh.add_attribute("vertex_normals")
        mesh.add_attribute("vertex_normals_shape")
        mesh.set_attribute("vertex_normals", vertex_normals)
        mesh.set_attribute("vertex_normals_shape", np.array(vertex_normals.shape))

    # get edge normals
    try:
        edge_normals = np.reshape(mesh.get_attribute("edge_normals"),
                                  np.int32(mesh.get_attribute("edge_normals_shape")))
        edge_map = np.reshape(np.int32(mesh.get_attribute("edge_map")),
                              np.int32(mesh.get_attribute("edge_map_shape")))
    except RuntimeError:
        edge_normals, _, edge_map = PyMesh.edge_normals(mesh.vertices, mesh.faces, face_normals)
        mesh.add_attribute("edge_normals")
        mesh.add_attribute("edge_map")
        mesh.add_attribute("edge_normals_shape")
        mesh.add_attribute("edge_map_shape")
        mesh.set_attribute("edge_normals", edge_normals)
        mesh.set_attribute("edge_map", edge_map)
        mesh.set_attribute("edge_normals_shape", np.array(edge_normals.shape))
        mesh.set_attribute("edge_map_shape", np.array(edge_map.shape))

    signed_distances, face_indices, closest_points, face_normals = bvh.lookup_signed(pts, face_normals, vertex_normals, edge_normals, edge_map);
    return signed_distances, face_indices, closest_points, face_normals;