Processing filters
The processing module of geos-mesh package contains filters to process meshes.
geos.mesh.processing.CreateConstantAttributePerRegion filter
CreateConstantAttributePerRegion is a vtk filter that allows to create an attribute with constant values per components for each chosen indexes of a reference/region attribute. If other region indexes exist values are set to nan for float type, -1 for int type or 0 for uint type.
Input mesh is either vtkMultiBlockDataSet or vtkDataSet and the region attribute must have one component. The relation index/values is given by a dictionary. Its keys are the indexes and its items are the list of values for each component. To use a handler of yours, set the variable ‘speHandler’ to True and add it using the member function addLoggerHandler.
By default, the value type is set to float32, their is one component and no name and the logger use an intern handler.
To use it:
from geos.mesh.processing.CreateConstantAttributePerRegion import CreateConstantAttributePerRegion
# Filter inputs.
mesh: Union[vtkMultiBlockDataSet, vtkDataSet]
regionName: str
dictRegionValues: dict[ Any, Any ]
newAttributeName: str
# Optional inputs.
valueNpType: type
nbComponents: int
componentNames: tuple[ str, ... ]
speHandler: bool
# Instantiate the filter
filter: CreateConstantAttributePerRegion = CreateConstantAttributePerRegion( mesh,
regionName,
dictRegionValues,
newAttributeName,
valueNpType,
nbComponents,
componentNames,
speHandler,
)
# Set your handler (only if speHandler is True).
yourHandler: logging.Handler
filter.addLoggerHandler( yourHandler )
# Do calculations.
filter.applyFilter()
- class geos.mesh.processing.CreateConstantAttributePerRegion.CreateConstantAttributePerRegion(mesh, regionName, dictRegionValues, newAttributeName, valueNpType=<class 'numpy.float32'>, nbComponents=1, componentNames=(), speHandler=False)[source]
Bases:
object
Create an attribute with constant value per region.
- Parameters:
mesh (Union[ vtkDataSet, vtkMultiBlockDataSet ]) – The mesh where to create the constant attribute per region.
regionName (str) – The name of the attribute with the region indexes.
dictRegionValues (dict[ Any, Any ]) – The dictionary with the region indexes as keys and the list of values as items.
newAttributeName (str) – The name of the new attribute to create.
valueNpType (type, optional) – The numpy scalar type for values. Defaults to numpy.float32.
nbComponents (int, optional) – Number of components for the new attribute. Defaults to 1.
componentNames (tuple[str,...], optional) – Name of the components for vectorial attributes. If one component, gives an empty tuple. Defaults to an empty tuple.
speHandler (bool, optional) – True to use a specific handler, False to use the internal handler. Defaults to False.
geos.mesh.processing.FillPartialArrays filter
Fill partial attributes of the input mesh with constant values per component.
Input mesh is vtkMultiBlockDataSet and attributes to fill must be partial.
The list of filling values per attribute is given by a dictionary. Its keys are the attribute names and its items are the list of filling values for each component.
If the list of filling value is None, attributes are filled with the same constant value for each component; 0 for uint data, -1 for int data and nan for float data.
To use a handler of yours for the logger, set the variable ‘speHandler’ to True and add it to the filter with the member function setLoggerHandler.
To use it:
from geos.mesh.processing.FillPartialArrays import FillPartialArrays
# Filter inputs.
multiBlockDataSet: vtkMultiBlockDataSet
dictAttributesValues: dict[ str, Union[ list[ Any ], None ] ]
# Optional inputs.
speHandler: bool
# Instantiate the filter.
filter: FillPartialArrays = FillPartialArrays( multiBlockDataSet, dictAttributesValues, speHandler )
# Set the handler of yours (only if speHandler is True).
yourHandler: logging.Handler
filter.setLoggerHandler( yourHandler )
# Do calculations.
filter.applyFilter()
- class geos.mesh.processing.FillPartialArrays.FillPartialArrays(multiBlockDataSet, dictAttributesValues, speHandler=False)[source]
Bases:
object
Fill partial attributes with constant value per component.
- If the list of filling values for an attribute is None, it will filled with the default value for each component:
0 for uint data. -1 for int data. nan for float data.
- Parameters:
multiBlockDataSet (vtkMultiBlockDataSet) – The mesh where to fill the attribute.
dictAttributesValues (dict[str, Any]) – The dictionary with the attribute to fill as keys and the list of filling values as items.
speHandler (bool, optional) – True to use a specific handler, False to use the internal handler. Defaults to False.
geos.mesh.processing.meshQualityMetricHelpers module
Helpers for MeshQuality metrics.
- geos.mesh.processing.meshQualityMetricHelpers.CELL_QUALITY_METRIC_ADDITIONAL_START_INDEX = 100
start index of additional cell quality metrics
- class geos.mesh.processing.meshQualityMetricHelpers.CellQualityMetricAdditionalEnum(value)[source]
Bases:
CellQualityMetricEnum
Additional cell quality metric enumeration.
The order of boolean is the same as getAllCellTypes method: (VTK_TRIANGLE, VTK_QUAD, VTK_TETRA, VTK_PYRAMID, VTK_WEDGE, VTK_HEXAHEDRON)
Metric index starts at 100 to prevent from conflicts with basic metrics and is incremented in the order of the enumeration.
Define the enumeration to add attributes to mesh quality measures.
- Parameters:
metricIndex (int) – Index of QualityMeasureTypes
name (str) – Name of the metric
applicableToCellTypes (tuple[bool, ...]) – Tuple defining for each cell type if the metric is applicable.
qualityRanges (tuple[QualityRange | None,...]) – Quality range limits for each cell type starting from best to worst quality.
- MAXIMUM_ASPECT_RATIO = (100, 'Maximum Aspect Ratio', (False, False, False, True, True, True), (None, None, None, QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf))))
maximum of aspect ratio over all tets .
- class geos.mesh.processing.meshQualityMetricHelpers.CellQualityMetricEnum(value)[source]
Bases:
MeshQualityMetricEnum
An enumeration.
Define the enumeration to add attributes to mesh quality measures.
- Parameters:
metricIndex (int) – Index of QualityMeasureTypes
name (str) – Name of the metric
applicableToCellTypes (tuple[bool, ...]) – Tuple defining for each cell type if the metric is applicable.
qualityRanges (tuple[QualityRange | None,...]) – Quality range limits for each cell type starting from best to worst quality.
- getApplicableCellTypes()[source]
Get the list of cell type indexes the metric applies to.
- Returns:
Set of cell type indexes
- Return type:
set[int]
- class geos.mesh.processing.meshQualityMetricHelpers.MeshQualityMetricEnum(value)[source]
Bases:
Enum
An enumeration.
Define the enumeration to add attributes to mesh quality measures.
- Parameters:
metricIndex (int) – Index of QualityMeasureTypes
name (str) – Name of the metric
- geos.mesh.processing.meshQualityMetricHelpers.QUALITY_METRIC_OTHER_START_INDEX = 200
start index of other mesh quality metrics
- class geos.mesh.processing.meshQualityMetricHelpers.QualityMetricOtherEnum(value)[source]
Bases:
MeshQualityMetricEnum
Additional metrics that apply to the mesh, not to specific cell type.
Metric index starts at 200 to prevent from conflicts with other metrics and is incremented in the order of the enumeration.
Define the enumeration to add attributes to mesh quality measures.
- Parameters:
metricIndex (int) – Index of QualityMeasureTypes
name (str) – Name of the metric
- INCIDENT_VERTEX_COUNT = (200, 'Incident Vertex Count')
number of incident edges for each vertex
- class geos.mesh.processing.meshQualityMetricHelpers.QualityRange(acceptableRange, normalRange, fullRange)[source]
Bases:
object
Defines metric quality ranges.
- acceptableRange
- fullRange
- normalRange
- class geos.mesh.processing.meshQualityMetricHelpers.VtkCellQualityMetricEnum(value)[source]
Bases:
CellQualityMetricEnum
Cell quality metric enumeration.
The order of boolean is the same as getAllCellTypes method: (VTK_TRIANGLE, VTK_QUAD, VTK_TETRA, VTK_PYRAMID, VTK_WEDGE, VTK_HEXAHEDRON)
Caution
The order of the enum must follow the one of vtkMeshQuality.QualityMeasureTypes.
Define the enumeration to add attributes to mesh quality measures.
- Parameters:
metricIndex (int) – Index of QualityMeasureTypes
name (str) – Name of the metric
applicableToCellTypes (tuple[bool, ...]) – Tuple defining for each cell type if the metric is applicable.
qualityRanges (tuple[QualityRange | None,...]) – Quality range limits for each cell type starting from best to worst quality.
- AREA = (28, 'Area (m2)', (True, True, False, False, False, False), (QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(-inf, inf)), None, None, None, None))
polygon area
- ASPECT_FROBENIUS = (3, 'Aspect Frobenius', (True, False, True, False, False, False), (QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), None, QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), None, None, None))
sum of the edge lengths squared divided by the area for triangles. Adapted for Tetrahedron. normalized so that equal to 1 when the element is equilateral triangle
- ASPECT_GAMMA = (27, 'Aspect Gamma', (False, False, True, False, False, False), (None, None, QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(0.0, inf)), None, None, None))
ratio of root-mean-square edge length to volume. normalizing the metric to a value of 1 for equilateral tetrahedra
- ASPECT_RATIO = (1, 'Aspect Ratio', (True, True, True, False, False, False), (QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf)), None, None, None))
ratio of the maximum edge length to the inradius (for simplicial elements but adapted for quads). may be adapted for polyhedron other than tet by splitting the polyhedron in tet and computing the max of aspect ratio of all tet
- COLLAPSE_RATIO = (7, 'Collapse Ratio', (False, False, True, False, False, False), (None, None, QualityRange(acceptableRange=(0.1, 1.0), normalRange=(0.0, inf), fullRange=(0.0, inf)), None, None, None))
the smallest ratio of the height of a vertex above its opposing triangle to the longest edge of that opposing triangle across all vertices of the tetrahedron
- CONDITION = (9, 'Condition', (True, True, True, False, True, True), (QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 4.0), normalRange=(1.0, 12.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf)), None, QualityRange(acceptableRange=(1.0, 4.0), normalRange=(1.0, 12.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 4.0), normalRange=(1.0, 12.0), fullRange=(1.0, inf))))
condition number of the weighted Jacobian matrix.
- DIAGONAL = (21, 'Diagonal', (False, False, False, False, False, True), (None, None, None, None, None, QualityRange(acceptableRange=(0.65, 1.0), normalRange=(0.0, 1.0), fullRange=(0.0, inf))))
ratio of the minimum diagonal length to the maximum diagonal length
- DIMENSION = (22, 'Dimension (m)', (False, False, False, False, False, True), (None, None, None, None, None, QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(0.0, inf))))
- DISTORTION = (15, 'Distortion', (True, True, True, False, True, True), (QualityRange(acceptableRange=(0.5, 1.0), normalRange=(0.0, 1.0), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.5, 1.0), normalRange=(0.0, 1.0), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.5, 1.0), normalRange=(0.0, 1.0), fullRange=(-inf, inf)), None, QualityRange(acceptableRange=(0.5, 1.0), normalRange=(0.0, 1.0), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.5, 1.0), normalRange=(0.0, 1.0), fullRange=(-inf, inf))))
measure of how well-behaved the mapping from parameter space to world coordinates is. ratio of the minimum of Jacobian determinant to cell area/volume
- EDGE_RATIO = (0, 'Edge Ratio', (True, True, True, False, True, True), (QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf)), None, QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf))))
ratio of cell longest and shortest edge lengths
- EQUIANGLE_SKEW = (29, 'Equiangle Skew', (True, True, True, True, True, True), (QualityRange(acceptableRange=(0.0, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.0, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.0, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.0, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.0, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.0, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0))))
maximum of ratio of angular deviation from ideal element
- EQUIVOLUME_SKEW = (30, 'Equivolume Skew', (False, False, True, False, False, False), (None, None, QualityRange(acceptableRange=(0.0, 0.3), normalRange=(0.0, 0.9), fullRange=(0.0, 1.0)), None, None, None))
maximum of ratio of volume deviation from ideal element
- JACOBIAN = (25, 'Jacobian', (False, True, True, True, True, True), (None, QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(-inf, inf))))
minimum determinant of the Jacobian matrix evaluated at each corner and the center of the element
- MAXIMUM_ANGLE = (8, 'Maximum Angle (°)', (True, True, False, False, False, False), (QualityRange(acceptableRange=(60.0, 90.0), normalRange=(60.0, 180.0), fullRange=(0.0, 180.0)), QualityRange(acceptableRange=(90.0, 135.0), normalRange=(90.0, 360.0), fullRange=(0.0, 360.0)), None, None, None, None))
maximum angle between two neighboring edges for polygons / faces for tetrahedron.
- MAXIMUM_ASPECT_FROBENIUS = (5, 'Maximum Aspect Frobenius', (False, True, False, False, True, True), (None, QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), None, None, QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf))))
maximum of Aspect Frobenius over all triangles of Quads / tetrahedra of hexahedron
- MAXIMUM_EDGE_RATIO = (16, 'Maximum Edge Ratio', (False, True, False, False, False, True), (None, QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), None, None, None, QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf))))
maximum of edge ratio over all triangles of the cell
- MAXIMUM_STRETCH = (31, 'Maximum Stretch', (False, False, False, False, True, False), (None, None, None, None, QualityRange(acceptableRange=(0.25, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, inf)), None))
maximum stretch over tetrahedra
- MEAN_ASPECT_FROBENIUS = (32, 'Mean Aspect Frobenius', (False, False, False, False, True, False), (None, None, None, None, QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf)), None))
mean of Aspect Frobenius over all triangles of Quads / tetrahedra of hexahedron
- MEAN_RATIO = (33, 'Mean Ratio', (False, False, True, False, False, False), (None, None, QualityRange(acceptableRange=(0.0, 0.3), normalRange=(0.0, 0.9), fullRange=(0.0, 1.0)), None, None, None))
ratio of tetrahedron volume over the volume of an equilateral tetrahedron with the same root mean squared edge length
- MEDIAN_ASPECT_FROBENIUS = (4, 'Med Aspect Frobenius', (False, True, False, False, False, True), (None, QualityRange(acceptableRange=(1.0, 1.3), normalRange=(1.0, 3.0), fullRange=(1.0, inf)), None, None, None, QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 3.0), fullRange=(9.0, inf))))
median of Aspect Frobenius over all triangles of quads / tetrahedra of hexahedron
- MINIMUM_ANGLE = (6, 'Minimum Angle (°)', (True, True, True, False, False, False), (QualityRange(acceptableRange=(30.0, 60.0), normalRange=(0.0, 60.0), fullRange=(0.0, 360.0)), QualityRange(acceptableRange=(45.0, 90.0), normalRange=(0.0, 90.0), fullRange=(0.0, 360.0)), QualityRange(acceptableRange=(40.0, np.float64(70.52877936550931)), normalRange=(0.0, np.float64(70.52877936550931)), fullRange=(0.0, 360.0)), None, None, None))
minimum angle between two neighboring edges for polygons / faces for tetrahedron.
- NODAL_JACOBIAN_RATIO = (34, 'Nodal Jacobian Ratio', (False, False, False, False, False, True), (None, None, None, None, None, QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(0.0, inf))))
ratio between the largest and smallest Jacobian determinant value
- NONE = (37, 'None', (False, False, False, False, False, False), (None, None, None, None, None, None))
no metric
- NORMALIZED_INRADIUS = (35, 'Normalized Inradius', (True, False, True, False, False, False), (QualityRange(acceptableRange=(0.15, 0.5), normalRange=(-1.0, 1.0), fullRange=(-1.0, 1.0)), None, QualityRange(acceptableRange=(0.15, 0.5), normalRange=(-1.0, 1.0), fullRange=(-1.0, 1.0)), None, None, None))
ratio of the minimum sub-triangle inner radius to the outer triangle radius
- ODDY = (23, 'Oddy', (False, True, False, False, False, True), (None, QualityRange(acceptableRange=(0.0, 0.5), normalRange=(0.0, 1.5), fullRange=(0.0, inf)), None, None, None, QualityRange(acceptableRange=(0.0, 0.5), normalRange=(0.0, 1.5), fullRange=(0.0, inf))))
measures the maximum deviation of the metric tensor at the corners of the quadrilateral. Maximum of oddy for hexahedron.
- RADIUS_RATIO = (2, 'Radius Ratio', (True, True, True, False, False, False), (QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf)), QualityRange(acceptableRange=(1.0, 3.0), normalRange=(1.0, 9.0), fullRange=(1.0, inf)), None, None, None))
ratio between the radius of the inscribed circle/sphere to the radius of the circum-circle/sphere normalized so that the ratio yields 1 for equilateral cell
- RELATIVE_SIZE_SQUARED = (12, 'Relative Size Squared', (True, True, True, False, False, True), (QualityRange(acceptableRange=(0.25, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.3, 0.6), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.3, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), None, None, QualityRange(acceptableRange=(0.5, 1.0), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0))))
the minimum of the ratio of cell area/volume to the average area/volume of an ensemble of cells and its inverse.
- SCALED_JACOBIAN = (10, 'Scaled Jacobian', (True, True, True, True, True, True), (QualityRange(acceptableRange=(0.5, np.float64(1.1547005383792515)), normalRange=(np.float64(-1.1547005383792515), np.float64(1.1547005383792515)), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.3, 1.0), normalRange=(-1.0, 1.0), fullRange=(-1.0, inf)), QualityRange(acceptableRange=(0.5, np.float64(0.7071067811865476)), normalRange=(np.float64(-0.7071067811865476), np.float64(0.7071067811865476)), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.5, 1.0), normalRange=(-1.0, 1.0), fullRange=(-1.0, inf)), QualityRange(acceptableRange=(0.5, 1.0), normalRange=(-1.0, 1.0), fullRange=(-1.0, inf)), QualityRange(acceptableRange=(0.5, 1.0), normalRange=(-1.0, 1.0), fullRange=(-1.0, inf))))
Jacobian divided by the product of the lengths of the longest edges normalized so that a unit equilateral triangle has value 1.
- SHAPE = (13, 'Shape', (True, True, True, True, True, True), (QualityRange(acceptableRange=(0.25, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.3, 0.6), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.3, 0.6), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.3, 0.6), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.3, 0.6), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.3, 0.6), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0))))
inverse of Condition (Polygons) / Jacobian ratio
- SHAPE_AND_SIZE = (14, 'Shape And Size', (True, True, True, False, False, True), (QualityRange(acceptableRange=(0.25, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.2, 0.4), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), QualityRange(acceptableRange=(0.2, 0.4), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), None, None, QualityRange(acceptableRange=(0.2, 0.4), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0))))
relative size squared times shape
- SHEAR = (11, 'Shear', (False, True, False, False, False, True), (None, QualityRange(acceptableRange=(0.3, 0.6), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), None, None, None, QualityRange(acceptableRange=(0.3, 0.6), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0))))
same as Scaled Jacobian
- SHEAR_AND_SIZE = (24, 'Shear And Size', (False, True, False, False, False, True), (None, QualityRange(acceptableRange=(0.2, 0.4), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), None, None, None, QualityRange(acceptableRange=(0.2, 0.4), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0))))
relative size squared times shear
- SKEW = (17, 'Skew', (False, True, False, False, False, True), (None, QualityRange(acceptableRange=(0.5, 1.0), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0)), None, None, None, QualityRange(acceptableRange=(0.0, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, 1.0))))
measures the angle (absolute value of the cosine) between the principal axes.
- SQUISH_INDEX = (36, 'Squish Index', (False, False, True, True, True, True), (None, None, QualityRange(acceptableRange=(0.0, 0.3), normalRange=(0.0, 0.9), fullRange=(0.0, 1.0)), None, None, None))
measure used to quantify how far a cell deviates from orthogonality with respect to its face maximum of sinus of the angle between the vector from polyhedron center and face center and face normal yields 0 if vectors are parallel, 1 if they are orthogonal
- STRETCH = (20, 'Stretch', (False, True, False, False, False, True), (None, QualityRange(acceptableRange=(0.25, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, inf)), None, None, None, QualityRange(acceptableRange=(0.25, 0.5), normalRange=(0.0, 1.0), fullRange=(0.0, inf))))
ratio of minimum edge length to longest diagonal length
- TAPER = (18, 'Taper', (False, True, False, False, False, True), (None, QualityRange(acceptableRange=(0.0, 0.7), normalRange=(0.0, 2.0), fullRange=(0.0, inf)), None, None, None, QualityRange(acceptableRange=(0.0, 0.5), normalRange=(0.0, 1.5), fullRange=(0.0, inf))))
maximum ratio of cross derivative magnitude to principal axis magnitude
- VOLUME = (19, 'Volume (m3)', (False, False, True, True, True, True), (None, None, QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(-inf, inf)), QualityRange(acceptableRange=(0.0, inf), normalRange=(0.0, inf), fullRange=(-inf, inf))))
polyhedron volume
- WARPAGE = (26, 'Warpage', (False, True, False, False, False, False), (None, QualityRange(acceptableRange=(0.0, 0.7), normalRange=(0.0, 2.0), fullRange=(0.0, inf)), None, None, None, None))
the cosine of the minimum dihedral angle formed by planes intersecting in diagonals (to the fourth power)
- geos.mesh.processing.meshQualityMetricHelpers.getAllCellTypes()[source]
Get all cell type ids.
- Returns:
Tuple containing cell type ids.
- Return type:
tuple[int,…]
- geos.mesh.processing.meshQualityMetricHelpers.getAllCellTypesExtended()[source]
Get all cell type ids.
- Returns:
Tuple containing cell type ids.
- Return type:
tuple[int,…]
- geos.mesh.processing.meshQualityMetricHelpers.getCellQualityMeasureFromCellType(cellType)[source]
Get the indexes of mesh quality metrics defined for triangles.
- Returns:
Set of possible indexes.
- Return type:
set[int]
- geos.mesh.processing.meshQualityMetricHelpers.getChildrenCellTypes(parent)[source]
Get children cell type ids from parent id.
- Returns:
Tuple containing cell type ids.
- Return type:
tuple[int,…]
- geos.mesh.processing.meshQualityMetricHelpers.getCommonPolygonQualityMeasure()[source]
Get the indexes of mesh quality metrics defined for both triangles and quads.
- Returns:
Set of possible indexes.
- Return type:
set[int]
- geos.mesh.processing.meshQualityMetricHelpers.getCommonPolyhedraQualityMeasure()[source]
Get the indexes of mesh quality metrics defined for tetrahedra, pyramids, wedges and hexahedra.
- Returns:
Set of possible indexes.
- Return type:
set[int]
- geos.mesh.processing.meshQualityMetricHelpers.getHexQualityMeasure()[source]
Get the indexes of mesh quality metrics defined for hexahedra.
- Returns:
Set of possible indexes.
- Return type:
set[int]
- geos.mesh.processing.meshQualityMetricHelpers.getPolygonCellTypes()[source]
Get polygonal cell type ids.
- Returns:
Tuple containing cell type ids.
- Return type:
tuple[int,…]
- geos.mesh.processing.meshQualityMetricHelpers.getPolyhedronCellTypes()[source]
Get polyhedra cell type ids.
- Returns:
Tuple containing cell type ids.
- Return type:
tuple[int,…]
- geos.mesh.processing.meshQualityMetricHelpers.getPyramidQualityMeasure()[source]
Get the indexes of mesh quality metrics defined for pyramids.
- Returns:
Set of possible indexes.
- Return type:
set[int]
- geos.mesh.processing.meshQualityMetricHelpers.getQuadQualityMeasure()[source]
Get the indexes of mesh quality metrics defined for quads.
- Returns:
Set of possible indexes.
- Return type:
set[int]
- geos.mesh.processing.meshQualityMetricHelpers.getQualityMeasureIndexFromName(name)[source]
Get quality metric index from name.
- Parameters:
name (str) – Name of quality measure
- Returns:
Index of quality measure
- Return type:
int
- geos.mesh.processing.meshQualityMetricHelpers.getQualityMeasureNameFromIndex(metricIndex)[source]
Get quality metric name from index.
- Parameters:
metricIndex (int) – Index of quality measure
- Returns:
Name of quality measure. Returns None if metricIndex is undefined.
- Return type:
str | None
- geos.mesh.processing.meshQualityMetricHelpers.getQualityMetricFromIndex(metricIndex)[source]
Get quality metric from its index.
- Parameters:
metricIndex (int) – Metric index
- Raises:
IndexError – Metric index is out of range
- Returns:
Quality metric
- Return type:
MeshQualityMetricEnum | None
- geos.mesh.processing.meshQualityMetricHelpers.getQualityMetricsOther()[source]
Get the set of indexes of other mesh quality metric.
- Returns:
Other mesh quality metric indexes
- Return type:
set[int]
- geos.mesh.processing.meshQualityMetricHelpers.getTetQualityMeasure()[source]
Get the indexes of mesh quality metrics defined for tetrahedra.
- Returns:
Set of possible indexes.
- Return type:
set[int]
geos.mesh.processing.SplitMesh filter
SplitMesh module is a vtk filter that split cells of a mesh composed of Tetrahedra, pyramids, and hexahedra.
Filter input and output types are vtkUnstructuredGrid.
To use the filter:
from geos.mesh.processing.SplitMesh import SplitMesh
# filter inputs
input :vtkUnstructuredGrid
# instanciate the filter
filter :SplitMesh = SplitMesh()
# set input data object
filter.SetInputDataObject(input)
# do calculations
filter.Update()
# get output object
output :vtkUnstructuredGrid = filter.GetOutputDataObject(0)
- class geos.mesh.processing.SplitMesh.SplitMesh[source]
Bases:
VTKPythonAlgorithmBase
SplitMesh filter splits each cell using edge centers.
- FillInputPortInformation(port, info)[source]
Inherited from VTKPythonAlgorithmBase::RequestInformation.
- Parameters:
port (int) – input port
info (vtkInformationVector) – info
- Returns:
1 if calculation successfully ended, 0 otherwise.
- Return type:
int
- RequestData(request, inInfoVec, outInfoVec)[source]
Inherited from VTKPythonAlgorithmBase::RequestData.
- Parameters:
request (vtkInformation) – request
inInfoVec (list[vtkInformationVector]) – input objects
outInfoVec (vtkInformationVector) – output objects
- Returns:
1 if calculation successfully ended, 0 otherwise.
- Return type:
int
- RequestDataObject(request, inInfoVec, outInfoVec)[source]
Inherited from VTKPythonAlgorithmBase::RequestDataObject.
- Parameters:
request (vtkInformation) – request
inInfoVec (list[vtkInformationVector]) – input objects
outInfoVec (vtkInformationVector) – output objects
- Returns:
1 if calculation successfully ended, 0 otherwise.
- Return type:
int