Source code for Muscat.Containers.Filters.FilterObjects

# -*- coding: utf-8 -*-
#
# This file is subject to the terms and conditions defined in
# file 'LICENSE.txt', which is part of this source code package.
#
from __future__ import annotations
from typing import Optional, List, Any, Union, Collection, Iterator, Callable, Tuple, Type, Sequence, Set
from functools import reduce
import numpy as np

from Muscat.Types import ArrayLike, MuscatIndex
from Muscat.Helpers.Decorators import froze_it
from Muscat.Helpers.Logger import Info, Warning
from Muscat.Containers.Mesh import Mesh
from Muscat.Containers.ElementsContainers import ElementsContainer
import Muscat.Containers.ElementsDescription as ED
from Muscat.Containers.Filters.FilterBase import FilterBase, ElementsSelection, Intersect1D, ZoneType, UnionOrdered1D

ZT_CENTER = "CENTER"
ZT_ALL_NODES = "ALL_NODES"
ZT_AT_LEAST_ONE = "AT_LEAST_ONE"


NT_ALL_NODES = "ALL_NODES"
NT_AT_LEAST_ONE = "AT_LEAST_ONE"

NM_ALL_NODES = "ALL_NODES"
NM_AT_LEAST_ONE = "AT_LEAST_ONE"

NF_ALL_NODES = "ALL_NODES"
NF_AT_LEAST_ONE = "AT_LEAST_ONE"


[docs]@froze_it class ElementFilter(FilterBase): """Class for element filtering by dimensionality, zone, mask, elementType, tag, external element filter and/or nodal filter Parameters ---------- dimensionality : list[int], optional the dimensionality of the element to be included in this filter, by default None possible option are: 0 1 2 3 elementType :Optional[Union[List[ED.ElementType], ED.ElementType]], optional the name of a element type or a list of element types to be included in this filter, by default None CENTER = if the the center of the element is inside the zone ALL_NODES = if all nodes of the element are inside the selection (zone, tag, mask...) AT_LEAST_ONE = if at least one node of the element is inside the selection (zone, tag, mask...) zonesTreatment : str, optional ["CENTER" | "ALL_NODES" | "AT_LEAST_ONE" ], by default "CENTER" [ ZT_CENTER | ZT_ALL_NODES | ZT_AT_LEAST_ONE ], by default ZT_CENTER nTagsTreatment : str, optional ["ALL_NODES" | "AT_LEAST_ONE" ], by default "ALL_NODES" [NT_ALL_NODES | NT_AT_LEAST_ONE ], by default NT_ALL_NODES nMaskTreatment : str, optional ["ALL_NODES" | "AT_LEAST_ONE" ], by default "ALL_NODES" [NM_ALL_NODES | NM_AT_LEAST_ONE ], by default NM_ALL_NODES nFilterTreatment : str, optional ["ALL_NODES" | "AT_LEAST_ONE" ], by default "ALL_NODES" [NF_ALL_NODES | NF_AT_LEAST_ONE ], by default NF_ALL_NODES """ def __init__(self, dimensionality: Optional[Union[Sequence[MuscatIndex], MuscatIndex]] = None, elementType: Optional[Union[List[ED.ElementType], ED.ElementType]] = None, zonesTreatment: str = ZT_CENTER, nTagsTreatment: str = NT_ALL_NODES, nMaskTreatment: str = NM_ALL_NODES, nFilterTreatment: str = NF_ALL_NODES, zone: Optional[Union[ZoneType, List[ZoneType]]] = None, eMask: Optional[ArrayLike] = None, nMask: Optional[ArrayLike] = None, eTag: Optional[Union[str, List[str]]] = None, nTag: Optional[Union[str, List[str]]] = None, eFilter: Optional[ElementFilter] = None, nFilter: Optional[FilterBase] = None ) -> None: super().__init__(zone=zone, eMask=eMask, nMask=nMask, eTag=eTag, nTag=nTag, eFilter=eFilter, nFilter=nFilter) self.eFilter:ElementFilter self.dimensionality = list() if dimensionality is not None: self.SetDimensionality(dimensionality) self.zonesTreatment = ZT_CENTER self.SetZonesTreatment(zonesTreatment) self.nTagsTreatment = NT_ALL_NODES self.SetNTagsTreatment(nTagsTreatment) self.nMaskTreatment = NM_ALL_NODES self.SetNMaskTreatment(nMaskTreatment) self.nFilterTreatment = NF_ALL_NODES self.SetNFilterTreatment(nFilterTreatment) self.elementTypes: Set[ED.ElementType] = set() if elementType is not None: self.SetElementType(elementType) self.withError = False self.zonesField = None # this is a hack
[docs] def IsEquivalent(self, other: Any) -> bool: """To check if 2 element filters are equivalent Parameters ---------- other : Any other object to check the equivalency Returns ------- bool True if the two filters are equal """ res = super().IsEquivalent(other) if not res: return False if id(self) == id(other): return True if len(self.dimensionality) != len(other.dimensionality): return False else: if np.any(np.sort(self.dimensionality) != np.sort(other.dimensionality)): return False if self.zonesTreatment != other.zonesTreatment: return False if self.elementTypes != other.elementTypes: return False return True
def __str__(self) -> str: """Return the string representation Returns ------- str Return the string representation """ res = "ElementFilter\n" res += f" dimensionality : {self.dimensionality}\n" res += f" elementTypes : {self.elementTypes}\n" res += f" eTags : {self.eTags}\n" res += f" nTags : {self.nTags}\n" res += f" nTagsTreatment : {self.nTagsTreatment}\n" res += f" eMask : {self.eMask}\n" res += f" nMask : {self.nMask}\n" res += f" nMaskTreatment : {self.nMaskTreatment}\n" res += f" zones : {self.zones}\n" res += f" zonesTreatment : {self.zonesTreatment}\n" res += f" eFilter : {self.eFilter}\n" res += f" nFilter : {self.nFilter}\n" res += f" nFilterTreatment : {self.nFilterTreatment}\n" return res
[docs] def SetZonesTreatment(self, zonesTreatment: str) -> None: """Set the way the elements are selected based on the position CENTER = if the the center of the element is inside the zone ALL_NODES = if all nodes of the element are inside the selection (zone, tag, mask...) AT_LEAST_ONE = if at least one node of the element is inside the selection (zone, tag, mask...) Parameters ---------- zonesTreatment : str ["CENTER" | "ALL_NODES" | "AT_LEAST_ONE" ], by default "CENTER" [ ZT_CENTER | ZT_ALL_NODES | ZT_AT_LEAST_ONE ], by default ZT_CENTER Raises ------ Exception if the string is not permitted """ if zonesTreatment in [ZT_CENTER, ZT_ALL_NODES, ZT_AT_LEAST_ONE]: self.zonesTreatment = zonesTreatment else: raise RuntimeError(f"Zone treatment not valid ({zonesTreatment}), possible options are : {ZT_CENTER}, {ZT_ALL_NODES}, {ZT_AT_LEAST_ONE}")
[docs] def SetDimensionality(self, dimensionality: Union[Sequence[MuscatIndex], MuscatIndex]) -> None: """ Set the dimensionality of the elements to be treated Parameters ---------- dimensionality : Union[int, List[int]] a int to be added to the list or a list of all the dimensionalities to take into account Raises ------ TypeError if the input is not a list or an int """ if dimensionality is not None: if isinstance(dimensionality, (list, np.ndarray, tuple)): self.dimensionality = [int(x) for x in dimensionality] elif isinstance(dimensionality, int): self.dimensionality.append(dimensionality) else: raise TypeError(f" dimensionality must be a list[int] or an int: is ({type(dimensionality)}): {dimensionality}") if len(self.dimensionality) > 0: if np.min(self.dimensionality) < 0: raise RuntimeError("dimensionality must be greater or equal to zero") if np.max(self.dimensionality) > 3: raise RuntimeError("dimensionality must be smaller than 4")
[docs] def SetElementType(self, elementTypes: Union[List[ED.ElementType], ED.ElementType]) -> None: """Set the names of the element types to be included in this filter Parameters ---------- elementTypes : Union[List[ED.ElementType], ED.ElementType] an element Type or the list of element types """ self.elementTypes = set() acceptedTypes = ED.ElementsInfo.keys() if isinstance(elementTypes, (ED.ElementType, str)): self.AddElementType(elementTypes) return for et in elementTypes: self.elementTypes.add(ED.ElementType(et))
[docs] def SetNTagsTreatment(self, nTagsTreatment: str) -> None: """Set the way the elements are selected based on the nTags ALL_NODES = if all nodes of the element are inside the selection (zone, tag, mask...) AT_LEAST_ONE = if at least one node of the element is inside the selection (zone, tag, mask...) Parameters ---------- nTagsTreatment : str ["ALL_NODES" | "AT_LEAST_ONE" ], by default "ALL_NODES" [NT_ALL_NODES | NT_AT_LEAST_ONE ], by default NT_ALL_NODES Raises ------ Exception if the nTagsTreatment string is not in ["ALL_NODES" | "AT_LEAST_ONE" ] """ if nTagsTreatment in [NT_ALL_NODES, NT_AT_LEAST_ONE]: self.nTagsTreatment = nTagsTreatment else: raise Exception(f"NTag treatment not valid ({nTagsTreatment}), possible options are : {NT_ALL_NODES}, {NT_AT_LEAST_ONE}")
[docs] def SetNMaskTreatment(self, nMaskTreatment: str) -> None: """Set the way the elements are selected based on the nodal mask nMask ALL_NODES = if all nodes of the element are inside the selection (zone, tag, mask...) AT_LEAST_ONE = if at least one node of the element is inside the selection (zone, tag, mask...) Parameters ---------- nMaskTreatment : str ["ALL_NODES" | "AT_LEAST_ONE"] Raises ------ Exception if the nMaskTreatment string is not in ["ALL_NODES", "AT_LEAST_ONE"] """ if nMaskTreatment in [NM_ALL_NODES, NM_AT_LEAST_ONE]: self.nMaskTreatment = nMaskTreatment else: raise Exception(f"nodal mask treatment not valid ({nMaskTreatment}), possible options are : ['{NM_ALL_NODES}', '{NM_AT_LEAST_ONE}']")
[docs] def SetNFilterTreatment(self, nFilterTreatment: str) -> None: """Set the way the elements are selected based on the nodal mask nMask ALL_NODES = if all nodes of the element are inside the selection (zone, tag, mask...) AT_LEAST_ONE = if at least one node of the element is inside the selection (zone, tag, mask...) Parameters ---------- nFilterTreatment : str ["ALL_NODES" | "AT_LEAST_ONE"] Raises ------ Exception if the nFilterTreatment string is not in ["ALL_NODES", "AT_LEAST_ONE"] """ if nFilterTreatment in [NF_ALL_NODES, NF_AT_LEAST_ONE]: self.nFilterTreatment = nFilterTreatment else: raise Exception(f"nodal mask treatment not valid ({nFilterTreatment}), possible options are : ['{NF_ALL_NODES}', '{NF_AT_LEAST_ONE}']")
[docs] def AddElementType(self, elementType: Union[ED.ElementType, str]) -> None: """Add an element type to be included in this filter Parameters ---------- elementType : str the name of an element type """ self.elementTypes.add(ED.ElementType(elementType))
def _CheckDimensionality_(self, elements: ElementsContainer) -> Union[bool, None]: """Internal function check if a type of element must be included based on the dimensionality criteria Parameters ---------- elements : ElementsContainer the incoming ElementsContainer Returns ------- Union[bool,None] True if this type of elements must be included False if this type of elements must be excluded None if this the filtering by dimensionality is not active """ if len(self.dimensionality) == 0: return None else: elementDimension = ED.dimensionality[elements.elementType] return elementDimension in self.dimensionality def _CheckElementTypes_(self, elements: ElementsContainer) -> Union[bool, None]: """Internal function check if a type of element must be included based on the elementType criteria Parameters ---------- elements : ElementsContainer the incoming ElementsContainer Returns ------- Union[bool,None] True if this type of elements must be included False if this type of elements must be excluded None if this the filtering by dimensionality is not active """ if len(self.elementTypes) == 0: return None if elements.elementType in self.elementTypes: return True else: return False def _CheckZones_(self, mesh: Mesh, elements: ElementsContainer) -> Union[np.ndarray, None]: """Internal function check if a type of element must be included based on the zone Parameters ---------- mesh : Mesh the incoming mesh elements : ElementsContainer the incoming ElementsContainer Returns ------- Union[np.array,None] np.ndarray of the indices inside the selection None if this the filtering by zones is not active """ if len(self.zones) == 0: return None numberOfObjects = elements.GetNumberOfElements() res = np.zeros(numberOfObjects, dtype=bool) if self.zonesTreatment == ZT_CENTER: from Muscat.Containers.MeshTools import GetElementsCenters centers = GetElementsCenters(nodes=mesh.nodes, elements=elements) for zone in self.zones: res2 = zone(centers) <= 0 np.logical_or(res, res2, out=res) else: for zone in self.zones: if self.zonesTreatment == ZT_ALL_NODES: z = zone(mesh.nodes) <= 0 res2 = np.sum(z[elements.connectivity], axis=1) == elements.GetNumberOfNodesPerElement() elif self.zonesTreatment == ZT_AT_LEAST_ONE: z = zone(mesh.nodes) <= 0 res2 = np.sum(z[elements.connectivity], axis=1) > 0 else: # pragma: no cover raise Exception(f"zonesTreatment unknown: ({self.zonesTreatment})") np.logical_or(res, res2, out=res) return np.where(res)[0].astype(MuscatIndex) def _CheckNTags_(self, mesh: Mesh, elements: ElementsContainer) -> Union[np.ndarray, None]: """Internal function check if a element must be included based on the nTags Parameters ---------- elements : ElementsContainer the incoming ElementsContainer Returns ------- Union[bool,None] np.ndarray of the indices inside the selection None if this the filtering by nTags is not active """ if len(self.nTags) == 0: return None nodalMask = np.zeros(mesh.GetNumberOfNodes(), dtype=bool) for tagName in self.nTags: if tagName in mesh.nodesTags: nodalMask[mesh.nodesTags[tagName].GetIds()] = True else:# pragma: no cover Info(f"Warning: Nodal Tag ('{tagName}') not present in the mesh") return self._CheckNMask_(elements, nodalMask, self.nTagsTreatment) def _CheckNMask_(self, elements: ElementsContainer, nodalMask: np.ndarray, treatment: str) -> np.ndarray: """Internal function check if a element must be included based on the nodal mask Parameters ---------- elements : ElementsContainer the incoming ElementsContainer nodalMask : np.ndarray a nodal mask (boolean array of size number of nodes) treatment : str [ "ALL_NODES", "AT_LEAST_ONE" ] Returns ------- np.ndarray np.ndarray of the indices inside the selection None if this the filtering by nTags is not active Raises ------ Exception if the treatment is not in [ "ALL_NODES", "AT_LEAST_ONE" ] """ if treatment == NT_ALL_NODES: res = np.sum(nodalMask[elements.connectivity], axis=1) == elements.GetNumberOfNodesPerElement() elif treatment == NT_AT_LEAST_ONE: res = np.sum(nodalMask[elements.connectivity], axis=1) > 0 else: # pragma: no cover raise Exception(f"treatment ({treatment}) unknown") return np.where(res)[0].astype(MuscatIndex, copy=False)
[docs] def GetIdsToTreat(self, mesh: Mesh, elementType: ED.ElementType) -> np.ndarray: """Get the entities selected by this filter Parameters ---------- elementType : str Elements type to treat Returns ------- Union np.ndarray The filtered entities """ elements = mesh.elements[elementType] if self._CheckDimensionality_(elements) == False: return np.array([], dtype=MuscatIndex) if self._CheckElementTypes_(elements) == False: return np.array([], dtype=MuscatIndex) # check all the tags exist elementTags = mesh.GetNamesOfElementTags() #check if the etag tags are present on the mesh tagsToTreat = [] for t in self.eTags: if t not in elementTags: Warning(f"Warning Element Tag ('{t}') not present in the mesh") continue tagsToTreat.append(t) # if only non existent tags are present the selection is empty if len(tagsToTreat) == 0 and len(tagsToTreat) != len(self.eTags): return np.array([], dtype=MuscatIndex) res = self._CheckTags_(elements.tags, tagsToTreat) if self._Done_(res): return res # pragma: no coverage res_nTags = self._CheckNTags_(mesh, elements) if res_nTags is not None and len(res_nTags) == 0: return res_nTags res = Intersect1D(res, res_nTags) if self._Done_(res): return res # pragma: no coverage res_zones = self._CheckZones_(mesh, elements) # this is a hack to be able to access the evaluated zone over the used points # for the moment this is not part of the public API self.zonesField = res_zones res = Intersect1D(res, res_zones) if self._Done_(res): return res # pragma: no coverage # Element Filter if self.eFilter is not None: resEF = self.eFilter.GetIdsToTreat(mesh, elementType) else: resEF = None res = Intersect1D(res, resEF) if self._Done_(res): return res # pragma: no coverage ### node filter ## if self.nFilter is not None: resNF = self.nFilter.GetNodesIndices(mesh) nMask = np.zeros(mesh.GetNumberOfNodes(), dtype=bool) nMask[resNF] = True res_nFilter = self._CheckNMask_(elements, nMask, self.nFilterTreatment) res = Intersect1D(res, res_nFilter) if self._Done_(res): return res # pragma: no coverage init = 0 for name, data in mesh.elements.items(): if name == elements.elementType: break init += data.GetNumberOfElements() res_mask = self._CheckMask_(self.eMask, init, elements.GetNumberOfElements()) res = Intersect1D(res, res_mask) if self._Done_(res): return res # pragma: no coverage if self.nMask is not None: res_nMask = self._CheckNMask_(elements, self.nMask, self.nMaskTreatment) res = Intersect1D(res, res_nMask) if self._Done_(res): return res # pragma: no coverage if res is None: return np.arange(elements.GetNumberOfElements(), dtype=MuscatIndex) else: return res
def __call__(self, mesh: Mesh) -> Iterator[ElementsSelection]: """Iterator over the selection :example: myFilter = ElementFilter(,dimensionality=2) for selection in myFilter(myMesh): print("This function is called only for 2D elements") print("Number of 2D element of type {selection.elements.elementType} is : {selection.Size()}" Parameters ---------- mesh : Mesh a mesh to filter Yields ------ Iterator[ElementsSelection] the selection for the current elementType Raises ------ Exception if withError is activated raise an exception to help the debugging """ meshOffset = 0 selectionOffset = 0 for elementType, elements in mesh.elements.items(): ids = self.GetIdsToTreat(mesh, elementType) if len(ids) > 0: yield ElementsSelection(mesh=mesh, indices=ids, elements=elements, meshOffset=meshOffset, selectionOffset=selectionOffset) meshOffset += elements.GetNumberOfElements() selectionOffset += len(ids) if selectionOffset == 0: if self.withError: raise Exception("Error!! Zero element in the element filter: \n" + str(self)) else: Info("Error!! Zero element in the element filter : \n" + str(self))
[docs] def GetNodesIndices(self, mesh: Mesh) -> np.ndarray: """Return the nodes selected by this filter Parameters ---------- mesh : Mesh the mesh Returns ------- np.ndarray the nodes indices selected by the filter """ res = np.empty((0), dtype=MuscatIndex) for selection in self(mesh): nodes = selection.elements.GetNodesIndexFor(selection.indices) res = UnionOrdered1D(res, nodes) return res
[docs]@froze_it class NodeFilter(FilterBase): """class for filtering nodes, bye , zone, mask, tag, external element filter and/or nodal filter the rest of the parameter are passed to the constructor of the base class Filter """ def __init__(self, zone: Optional[Union[ZoneType, List[ZoneType]]] = None, eMask: Optional[ArrayLike] = None, nMask: Optional[ArrayLike] = None, eTag: Optional[Union[str, List[str]]] = None, nTag: Optional[Union[str, List[str]]] = None, eFilter: Optional[FilterBase] = None, nFilter: Optional[FilterBase] = None ) -> None: super().__init__(zone=zone, eMask=eMask, nMask=nMask, eTag=eTag, nTag=nTag, eFilter=eFilter, nFilter=nFilter) def _CheckZones_(self, pos: ArrayLike) -> Union[np.ndarray, None]: """Internal function to compute the indices to be treated based on the zones Parameters ---------- pos : ArrayLike (n,3) size array with the positions to be treated numberOfObjects : MuscatIndex total number of points Returns ------- Union[np.ndarray,None] list of nodes to treat or None for all nodes """ if len(self.zones) == 0: return None if len(self.zones) == 1: return np.where(self.zones[0](pos) <= 0)[0].astype(MuscatIndex) mask = reduce(np.logical_or, (zone(pos) <= 0 for zone in self.zones)) return np.where(mask)[0].astype(MuscatIndex)
[docs] def GetNodesIndices(self, mesh: Mesh) -> np.ndarray: """Get the nodes indices selected by this filter Parameters ---------- mesh : Mesh mesh to apply the filter Returns ------- Union np.ndarray The filtered nodes indices """ if len(self.eTags) > 0: resET = ElementFilter(eTag=self.eTags).GetNodesIndices(mesh) else: resET = None if self._Done_(resET): return resET else: res = resET if self.eFilter is not None: resEF = self.eFilter.GetNodesIndices(mesh) else: resEF = None res = Intersect1D(res, resEF) if self._Done_(res): return res # pragma: no coverage if self.eMask is not None: resEM = ElementFilter(eMask=self.eMask).GetNodesIndices(mesh) else: resEM = None res = Intersect1D(res, resEM) if self._Done_(res): return res # pragma: no coverage # check all the tags exist nodalTagNames = mesh.nodesTags.keys() for t in self.nTags: if t not in nodalTagNames: raise RuntimeError(f"Nodal Tag ({t}) not present in the mesh, nodal tags in the mesh are {nodalTagNames}") resT = self._CheckTags_(mesh.nodesTags, self.nTags) res = Intersect1D(res, resT) if self._Done_(res): return res # pragma: no coverage resZ = self._CheckZones_(mesh.nodes) res = Intersect1D(res, resZ) if self._Done_(res): return res # pragma: no coverage resM = self._CheckMask_(self.nMask, 0, mesh.GetNumberOfNodes()) res = Intersect1D(res, resM) if self._Done_(res): return res # pragma: no coverage if res is None: return np.arange(mesh.GetNumberOfNodes(), dtype=MuscatIndex) else: return res
[docs]def CheckIntegrity(GUI: bool = False) -> str: from Muscat.Containers import MeshCreationTools as UMCT from Muscat.Helpers.CheckTools import MustFailFunction, MustFailFunctionWith from Muscat.Containers.Filters.FilterTools import GetFrozenFilter mesh = UMCT.CreateSquare(dimensions=[10, 10], origin=[0, 0], spacing=[10/9, 10/9], ofTriangles=True) mesh.nodesTags.CreateTag("first Points").SetIds([0]) mesh.nodesTags.CreateTag("Second Points").SetIds([1]) mesh.elements.GetElementsOfType(ED.Triangle_3).tags.CreateTag("first Elem").SetIds([0]) mesh.elements.GetElementsOfType(ED.Triangle_3).tags.CreateTag("second Elem").SetIds([1]) mesh.elements.GetElementsOfType(ED.Triangle_3).tags.CreateTag("empty") numberOfTriangles = mesh.elements.GetElementsOfType(ED.Triangle_3).GetNumberOfElements() ############### Node Filter ################# nf = NodeFilter() np.testing.assert_equal(nf.GetNodesIndices(mesh), np.arange(100)) nf.SetETags(["first Elem"]) indices = nf.GetNodesIndices(mesh) np.testing.assert_equal(indices, [0, 10, 11]) nf.SetETags(["empty"]) indices = nf.GetNodesIndices(mesh) np.testing.assert_equal(indices, []) nf.SetETags([]) nf.eFilter = ElementFilter(eTag=["first Elem"]) indices = nf.GetNodesIndices(mesh) np.testing.assert_equal(indices, [0, 10, 11]) nf.eFilter = None nf.eMask = np.zeros(mesh.elements.GetElementsOfType(ED.Triangle_3).GetNumberOfElements(), dtype=bool) nf.eMask[3] = True # this is a bar indices = nf.GetNodesIndices(mesh) np.testing.assert_equal(indices, [2, 3]) nf.eMask = None nf.SetZones([lambda xyz: xyz[:, 0] - 0.0001]) indices = nf.GetNodesIndices(mesh) np.testing.assert_equal(indices, [0, 1, 2, 3, 4, 5, 6, 7, 8, 9]) nf.SetZones([lambda xyz: xyz[:, 0] - 0.0001, lambda xyz: xyz[:, 1] - 0.0001, ]) indices = nf.GetNodesIndices(mesh) np.testing.assert_equal(indices, [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 20, 30, 40, 50, 60, 70, 80, 90]) MustFailFunction(GetFrozenFilter, NodeFilter(nTag="non existent tag"), mesh, onNodes=True, onElements=False) ############### Elements Filter ################# ef = ElementFilter(dimensionality=0, elementType=ED.Triangle_3) ef.AddElementType(ED.Triangle_3) ef.AddElementType(ED.ElementType["Triangle_3"]) MustFailFunction(ElementFilter, dimensionality=3.4) MustFailFunction(ElementFilter, elementType="mySuperType") MustFailFunctionWith(ValueError, ElementFilter, elementType=["mySuperType"]) MustFailFunctionWith(RuntimeError, ElementFilter, zonesTreatment="mySuperTreatment") MustFailFunctionWith(Exception, ElementFilter, nTagsTreatment="mySuperTreatment") ef = ElementFilter(dimensionality=[2], elementType=[ED.Triangle_3]) assert (ef.IsEquivalent(ef)) assert (ef.IsEquivalent(0) == False) assert (ef.IsEquivalent(ElementFilter()) == False) assert (ef.IsEquivalent(ElementFilter(dimensionality=0)) == False) assert (ef.IsEquivalent(ElementFilter(dimensionality=2, elementType=ED.Bar_2)) == False) ef.zones = [lambda xyz: xyz[:, 0] - 0.0001, lambda xyz: xyz[:, 1] - 0.0001] assert (ef.IsEquivalent(ElementFilter(dimensionality=2, zonesTreatment=ZT_ALL_NODES, zone=ef.zones)) == False) ef.zones = [] assert (ef.IsEquivalent(ElementFilter(dimensionality=[2], elementType=[ED.Triangle_3]))) for selection in ef(mesh): # assert(selection.indices) assert (selection.Size() == numberOfTriangles) for selection in ElementFilter(dimensionality=[2], elementType=[ED.Bar_2])(mesh): raise Exception("normally the filter must be empty ") # pragma: no cover np.testing.assert_equal(GetFrozenFilter(ElementFilter(zone=lambda p: (p[:, 0]-0.2), zonesTreatment=ZT_ALL_NODES), mesh).GetElementSelectionSize(), 9) np.testing.assert_equal(GetFrozenFilter(ElementFilter(zone=lambda p: (p[:, 0]-0.001), zonesTreatment=ZT_AT_LEAST_ONE), mesh).GetElementSelectionSize(), 29) np.testing.assert_equal(GetFrozenFilter(ElementFilter(zone=lambda p: (p[:, 0]-0.2), zonesTreatment=ZT_CENTER), mesh).GetElementSelectionSize(), 9) np.testing.assert_equal(GetFrozenFilter(ElementFilter(nTag="Second Points", nTagsTreatment=NT_ALL_NODES), mesh).GetElementSelectionSize(), 0) np.testing.assert_equal(GetFrozenFilter(ElementFilter(nTag="Second Points", nTagsTreatment=NT_AT_LEAST_ONE), mesh).GetElementSelectionSize(), 5) nMask = np.zeros(mesh.GetNumberOfNodes(), dtype=bool) nMask[0] = True np.testing.assert_equal(GetFrozenFilter(ElementFilter(nMask=nMask, nMaskTreatment=NM_AT_LEAST_ONE), mesh).GetElementSelectionSize(), 4) np.testing.assert_equal(GetFrozenFilter(ElementFilter(dimensionality=0), mesh).GetElementSelectionSize(), 0) np.testing.assert_equal(GetFrozenFilter(ElementFilter(elementType=[ED.Bar_2]), mesh).GetElementSelectionSize(), 36) mesh.GetElementsOfType(ED.Bar_2) # add an empty bar container # from Muscat.Actions.OpenInParaView import OpenInParaView # OpenInParaView(mesh,wait=True) # from Muscat.Bridges.vtkBridge import PlotMesh # PlotMesh(mesh) np.testing.assert_equal(GetFrozenFilter(ElementFilter(elementType=[ED.Triangle_3]), mesh).GetElementSelectionSize(), 162) np.testing.assert_equal(GetFrozenFilter(ElementFilter(elementType=[ED.Bar_2], nFilter=NodeFilter(nTag="x0y0"), nFilterTreatment=NF_AT_LEAST_ONE), mesh).GetElementSelectionSize(), 2) def TestEmpty(filter, mesh): for selection in filter(mesh): pass # pragma: no cover ef = ElementFilter(elementType=[ED.Bar_3]) ef.withError = True MustFailFunction(TestEmpty, ef, mesh) MustFailFunction(ef.SetNMaskTreatment, ZT_CENTER) MustFailFunction(ef.SetNFilterTreatment, ZT_CENTER) MustFailFunction(ef.SetDimensionality, "0") MustFailFunction(ef.SetDimensionality, 4) MustFailFunction(ef.SetDimensionality, -1) GetFrozenFilter(ElementFilter(eTag="non existent tag"), mesh) np.testing.assert_equal(GetFrozenFilter(ElementFilter(eFilter=ElementFilter()), mesh).GetElementSelectionSize(), 198) return "ok"
if __name__ == '__main__': print(CheckIntegrity(GUI=True)) # pragma: no cover