Source code for Muscat.Bridges.PyVistaBridge

#
# This file is subject to the terms and conditions defined in
# file 'LICENSE.txt', which is part of this source code package.
#

from typing import Union
from Muscat.MeshContainers.Mesh import Mesh

[docs] def MeshToPyVista(mesh: Mesh, tagsAsFields: bool=False): from Muscat.Bridges.vtkBridge import MeshToVtk import pyvista as pv return pv.wrap(MeshToVtk(mesh, tagsAsFields=tagsAsFields))
[docs] def PyVistaToMesh(pvmesh, fieldsAsTags:bool=False) -> Mesh: from Muscat.Bridges.vtkBridge import VtkToMesh return VtkToMesh(pvmesh, fieldsAsTags=fieldsAsTags)
[docs] def PlotMesh(mesh: Union[Mesh,'DataSet'], tagsAsFields: bool = False, scalars= None, show_edges=None, cpos= None, show_scalar_bar=None, **pyVista_Plotter_options): # pragma: no cover """Plot a mesh using PyVista, the user can press "n" to plot the next field and "c" to cycle thro the components in the field. Parameters ---------- mesh : Mesh, vtk mesh or a pyvista mesh the mesh to be plotted tagsAsFields : bool, optional if True the tags of the mesh (in the case of a Muscat Mesh) are converted to fields, by default False pyVista_Plotter_options : All the rest of the kwargs are passed to the pyvista Plotter (like theme and other options) Examples -------- Simply show the plot of a mesh. >>> from Muscat.MeshTools.MeshCreationTools import CreateDisk >>> from pyvista.plotting.themes import DarkTheme >>> from Muscat.Simple import PlotMesh Mesh creation, and fill some synthetic data >>> mesh=CreateDisk() >>> mesh.nodeFields = {'x': mesh.nodes[:, 0].flatten(), 'Pos': mesh.nodes} Plot the mesh >>> PlotMesh(n, theme=DarkTheme()) """ if isinstance(mesh, Mesh): pyVistaMesh = MeshToPyVista(mesh, tagsAsFields=tagsAsFields) else: pyVistaMesh = mesh import pyvista as pv pl = pv.Plotter(**pyVista_Plotter_options) current = 0 component = 0 nbpd = list(pyVistaMesh.point_data.keys()) nbcd = list(pyVistaMesh.cell_data.keys()) if scalars is None and len(nbpd): scalars=nbpd[0] actor = pl.add_mesh(pyVistaMesh, scalars=scalars, show_edges=show_edges, show_scalar_bar=show_scalar_bar ) class Handler: def __init__(self): self.component = 0 self.current = 0 def nextComponent(self): self.component += 1 pl.render() self.render_handler() def toggle(self): self.current += 1 self.render_handler() def OpenInParaView(self): if isinstance(mesh, Mesh): from Muscat.Actions.OpenInParaView import OpenInParaView OpenInParaView(mesh) def render_handler(self): i = self.current%len(nbpd+nbcd) if i < len(nbpd): preference = 'point' scalar_name = nbpd[i] scalar_value = pyVistaMesh.point_data[scalar_name] else: preference = 'cell' scalar_name = nbcd[i-len(nbpd)] scalar_value = pyVistaMesh.cell_data[scalar_name] if len(scalar_value.shape)>1: local_component = self.component%scalar_value.shape[1] scalar_name += f"_{local_component}" else: local_component = 0 scalar_name += f" (on {preference})" actor.mapper.set_scalars(scalar_value.astype(float)*0.5, scalars_name=scalar_name, preference=None, component=local_component) pl.scalar_bar.SetTitle(scalar_name) pl.render() keyHandler = Handler() pl.add_key_event('n', keyHandler.toggle) pl.add_key_event('c', keyHandler.nextComponent) pl.add_key_event('o', keyHandler.OpenInParaView) pl.add_text("n: Next field c:Next component, w:wireframe, s:solid, r: reset, o:OpenInParaView, q:quit",font_size=7, position='lower_left') pl.show(cpos= cpos)
[docs] def ReadMeshWithPyVista(filename, fieldsAsTags=False)-> Mesh: """ Function to delegate directly the reading of the file to Pyvista https://docs.pyvista.org/api/utilities/_autosummary/pyvista.read.html and then cast pyvista mesh to muscat mesh """ import pyvista as pv pvmesh = pv.read(filename) return PyVistaToMesh(pvmesh, fieldsAsTags=fieldsAsTags)
[docs] def CheckIntegrity(GUI: bool = False): from Muscat.Helpers.CheckTools import SkipTest if SkipTest("PYVISTA_NO_FAIL"): # pragma: no cover return "skip" try: import pyvista except: # pragma: no cover return "skip : pyvista not installed" import Muscat.MeshContainers.ElementsDescription as ED from Muscat.MeshTools.MeshCreationTools import CreateMeshOf from Muscat.MeshTools.MeshTools import IsClose points = [[0, 0, 0], [1, 0, 0], [0, 1, 0], [0, 0, 1]] tet = [[0, 1, 2, 3]] res = CreateMeshOf(points, tet, elemName=ED.Tetrahedron_4) res.nodeFields = {"x": res.nodes[:, 0].flatten(), "Pos": res.nodes} res.nodesTags.CreateTag("FirstPoint").AddToTag(0) res.elemFields = {"SecondPoint": res.GetElementsOfType(ED.Tetrahedron_4).connectivity[:, 1].flatten().astype(float), "conn": res.GetElementsOfType(ED.Tetrahedron_4).connectivity} res.GetElementsOfType(ED.Tetrahedron_4).tags.CreateTag("FirstTetrahedron").AddToTag(0) sol = MeshToPyVista(res, tagsAsFields=True) resII = PyVistaToMesh(sol, fieldsAsTags=True) print(res) print(resII) if not IsClose(res, resII): # pragma: no cover raise (Exception("The meshes are not equal")) if GUI: # pragma: no cover PlotMesh(res) PlotMesh(resII, eye_dome_lighting=True, cpos=[-1, -1, 0.2], color=True) return "ok"
if __name__ == "__main__": print(CheckIntegrity(True)) # pragma: no cover