Examples of usage of filtersΒΆ
Here we show how to use node and element filters. The main Muscat modules used in this example are listed below:
ElementFilter
: allows for filtering the elements within your mesh.NodeFilter
: allows for filtering the nodes within your mesh.IntersectionFilter
: useful for combining intersecting filters together.UnionFilter
: useful for combining filters together.
[1]:
from Muscat.Containers.Filters.FilterTools import ElementFilter, NodeFilter
from Muscat.Containers.Filters.FilterOperators import UnionFilter, IntersectionFilter
import numpy as np
from scipy.spatial import Delaunay
from Muscat.Containers import MeshCreationTools as UMCT
import pyvista
pyvista.global_theme._jupyter_backend = 'panel' # remove this line to get interactive 3D plots
def PlotnTag(mesh, tag):
plotter = pv.Plotter()
actor = plotter.add_mesh(MeshToPyVista(mesh,TagsAsFields=True),show_edges=True, edge_color="grey", style="wireframe")
actor = plotter.add_mesh(MeshToPyVista(mesh,TagsAsFields=True), scalars=tag,show_edges=True, edge_color="grey", style="points")
plotter.enable_parallel_projection()
plotter.show(cpos="xy")
def PloteTag(mesh, tag):
plotter = pv.Plotter()
actor = plotter.add_mesh(MeshToPyVista(mesh,TagsAsFields=True),scalars=tag, show_edges=True, edge_color="grey", style="surface")
#actor = plotter.add_mesh(MeshToPyVista(mesh,TagsAsFields=True), scalars=tag,show_edges=True, edge_color="grey", style="points")
plotter.enable_parallel_projection()
plotter.show(cpos="xy")
For this simple example, we generate a dummy mesh using scipy.spatial.Delaunay
, and create a Mesh
object using Muscat.Containers.MeshCreationTools
.
[2]:
def create_mesh_example():
"""This function creates a dummy mesh using Delaunay triangulation.
"""
# Create a grid of points
x, y = np.meshgrid(np.linspace(0,1,100), np.linspace(0,1,100))
points = np.stack([x.ravel(), y.ravel()], axis=1)
# Generate the triangles with Delaunay
tri = Delaunay(points)
triangles = tri.simplices
# Create a Muscat Mesh using the CreateMeshOfTriangles utility
mesh = UMCT.CreateMeshOfTriangles(points, triangles)
return mesh
mesh = create_mesh_example()
[3]:
from Muscat.Bridges.PyVistaBridge import MeshToPyVista
import pyvista as pv
PloteTag(mesh,"2D")
Node filters can be used to define sets of nodes satisfying an inequality constraint of the form :math:f(x) < 0
.
In the example below, we select nodes that lie within the ball centered at :math:(0.5, 0.5)
and with radius :math:0.25
.
The node filter can then be used to get the indices of the nodes satisfying the constraint, which in turn can be used to create a new nodal tag.
[4]:
nf = NodeFilter(zone=lambda p: (np.linalg.norm(p - 0.5, axis=1) - 0.25))
indices = nf.GetNodesIndices(mesh)
mesh.GetNodalTag("circle").AddToTag(indices)
PlotnTag(mesh, "circle")
Multiple filters can be combined together using UnionFilter
.
[5]:
nf1 = NodeFilter(zone=lambda p: (np.linalg.norm(p - 0.6, axis=1) - 0.25))
nf2 = NodeFilter(zone=lambda p: (np.linalg.norm(p - 0.4, axis=1) - 0.25))
nf = UnionFilter(filters=[nf1, nf2])
indices = nf.GetNodesIndices(mesh)
mesh.GetNodalTag("two_circles_union").AddToTag(indices)
PlotnTag(mesh, "two_circles_union")
We can also compute the intersection of multiple filters using IntersectionFilter
.
[6]:
nf1 = NodeFilter(zone=lambda p: (np.linalg.norm(p - 0.6, axis=1) - 0.25))
nf2 = NodeFilter(zone=lambda p: (np.linalg.norm(p - 0.4, axis=1) - 0.25))
nf = IntersectionFilter(filters=[nf1, nf2])
indices = nf.GetNodesIndices(mesh)
mesh.GetNodalTag("two_circles_intersection").AddToTag(indices)
PlotnTag(mesh, "two_circles_intersection")
The same logic applies to ElementFilter
in Muscat.
[7]:
ne = ElementFilter(zone=lambda p: (np.linalg.norm(p - 0.5, axis=1) - 0.25))
for selection in ne(mesh):
selection.elements.tags.CreateTag("circle_element_filter").SetIds(selection.indices)
#indices = ne.GetNodesIndices(mesh)
#mesh.GetNodalTag("circle_element_filter").AddToTag(indices)
PloteTag(mesh, "circle_element_filter")
Element filters can also be combined together in a similar fashion.
Elements and nodes filters can be used together to define new sets of nodes. As an illustration, applying the ComputeSkin
method to our mesh creates a new mesh with an additional element tag Skin
. The latter corresponds to the set of elements that lie on the boundary of the mesh.
[8]:
from Muscat.Containers import MeshModificationTools as MMT
new_mesh = MMT.ComputeSkin(mesh, inPlace=True)
#Printing the mesh shows all the tags that we added previsouly and the new `Skin` tag for elements.
print(new_mesh)
Mesh
Number Of Nodes : 10000
Tags : (circle:1928) (two_circles_union:3230) (two_circles_intersection:624)
Number Of Elements : 19998
ElementsContainer, Type : (ElementType.Bar_2,396), Tags : (Skin:396)
ElementsContainer, Type : (ElementType.Triangle_3,19602), Tags : (2D:19602) (circle_element_filter:3827)
We can now easily determine the nodes that lie on the top of the mesh using a NodeFilter
that only considers the elements belonging to Skin
.
[9]:
nf = NodeFilter(eTag="Skin", zone = lambda p: -p[:, 1] + 0.99)
indices = nf.GetNodesIndices(new_mesh)
new_mesh.GetNodalTag("top").AddToTag(indices)
The mesh now has a new nodal tag top
:
[10]:
print(new_mesh)
PlotnTag(new_mesh,'top')
Mesh
Number Of Nodes : 10000
Tags : (circle:1928) (two_circles_union:3230) (two_circles_intersection:624) (top:100)
Number Of Elements : 19998
ElementsContainer, Type : (ElementType.Bar_2,396), Tags : (Skin:396)
ElementsContainer, Type : (ElementType.Triangle_3,19602), Tags : (2D:19602) (circle_element_filter:3827)