Utilities classes for processing filters
The tools folder contains utilities classes that are used in some of the processing filters.
FaultGeometry
- class geos.processing.tools.FaultGeometry.FaultGeometry(mesh, faultValues, faultAttribute, volumeMesh, outputDir='.', logger=None)[source]
Bases:
objectHandles fault surface extraction and normal computation with optimizations.
Initialize fault geometry with pre-computed topology.
- Parameters:
mesh (vtkUnstructuredGrid) – Pre-simulation mesh
faultValues (list[int]) – Values of fault attribute to consider
faultAttribute (str) – Fault attribute name in the mesh
volumeMesh (vtkUnstructuredGrid) – Volumic mesh
outputDir (str, optional) – Output directory Defaults is “.”.
logger (Union[Logger, None], optional) – A logger to manage the output messages. Defaults to None, an internal logger is used.
- computeDipStrikeFromCellBase(normals, tangent1, tangent2)[source]
Compute dip and strike angles from cell normal and tangent vectors.
- Assumptions:
Coordinate system: X = East, Y = North, Z = Up.
Vectors are provided per cell (shape: (nCells, 3)).
Input vectors are assumed to be orthonormal (n = t1 × t2).
- Parameters:
normals (npt.NDArray[np.float64]) – Normal vectors
tangent1 (npt.NDArray[np.float64]) – First tangent vector
tangent2 (npt.NDArray[np.float64]) – Second tangent vector
- Returns:
Dip and strike.
- Return type:
tuple[npt.NDArray[np.float64], npt.NDArray[np.float64]]
- diagnoseNormals()[source]
Diagnostic visualization to check normal quality.
Shows orthogonality and orientation issues.
- Returns:
The input surface, possibly annotated with an additional field indicating orthogonality errors.
- Return type:
vtkUnstructuredGrid
- getContributingCells(side='all')[source]
Get the extracted contributing cells.
- Parameters:
side (str) – ‘all’, ‘plus’, or ‘minus’
- Returns:
Contributing volume cells
- Return type:
vtkUnstructuredGrid
- Raises:
ValueError – Contributed cells not computed
ValueError – Invalid requested side
- getGeometricProperties()[source]
Get pre-computed geometric properties.
- Properties include:
‘volumes’: ndarray of cell volumes
‘centers’: ndarray of cell centers (nCells, 3)
‘distances’: ndarray of distances to nearest fault cell
‘faultTree’: KDTree for fault surface
- Returns:
geometric properties
- Return type:
dict[ str, Any ]
- initialize(processFaultsSeparately=True, saveContributionCells=True)[source]
One-time initialization: compute normals, adjacency topology, and geometric properties.
- Parameters:
processFaultsSeparately (bool) – Flag to process faults separately or not. Defaults is True.
saveContributionCells (bool) – Save the contributing cells as VTU. Defaults is True.
FaultVisualizer
- class geos.processing.tools.FaultVisualizer.Visualizer(profileSearchRadius=None, minDepthProfiles=None, maxDepthProfiles=None, savePlots=True, logger=None)[source]
Bases:
objectVisualization utilities.
Visualization utilities.
- Parameters:
profileSearchRadius (float | None) – Searching radius for determination of profile.
minDepthProfiles (float | None) – Minimum depth profile.
maxDepthProfiles – (float | None): Maximum depth profile.
savePlots (bool) – Flag to save the figures. Defaults is True,
logger (Union[Logger, None], optional) – A logger to manage the output messages. Defaults to None, an internal logger is used.
- plotDepthProfiles(surface, time, faultAttribute, path=PosixPath('.'), save=True, profileStartPoints=None, maxProfilePoints=1000, referenceProfileId=1)[source]
Plot vertical profiles along the fault showing stress and SCU vs depth.
- Parameters:
surface (vtkUnstructuredGrid) – Fault mesh.
time (float) – Time.
faultAttribute (str) – Fault attribute name.
path (Path) – Saving path. Defaults is current directory.
save (bool) – Flag to save plots. Defaults is True.
profileStartPoints (list[ tuple[ float, float ] ] | None) – List of start points for profile analysis
maxProfilePoints (int) – Max profile points displayed. Defaults is 1000.
referenceProfileId (int) – Id of profile to plot. Defaults is 1.
- plotMohrCoulombDiagram(surface, time, path=PosixPath('.'), save=True)[source]
Create Mohr-Coulomb diagram with depth coloring.
- Parameters:
surface (vtkUnstructuredGrid) – Fault mesh containing mohr attributes.
time (float) – Time.
path (Path) – Saving path.
save (bool) – Flag to save figures. Defaults is True.
- plotVolumeStressProfiles(volumeMesh, faultSurface, time, path=PosixPath('.'), save=True, profileStartPoints=None, maxProfilePoints=1000)[source]
Plot stress profiles in volume cells adjacent to the fault.
Extracts profiles through contributing cells on BOTH sides of the fault Shows plus side and minus side on the same plots for comparison.
- Parameters:
volumeMesh (vtkUnstructuredGrid) – Volumic mesh.
faultSurface (vtkUnstucturedGrid) – Fault mesh.
time (float) – Time.
path (Path) – Saving path. Defaults is current directory.
save (bool) – Flag to save plots. Defaults is True.
profileStartPoints (list[ tuple[ float, float ] ] | None) – List of start points for profile analysis
maxProfilePoints (int) – Max profile points displayed. Defaults is 1000.
MohrCoulomb
- class geos.processing.tools.MohrCoulomb.MohrCoulombAnalysis(surface, cohesion, frictionAngle, logger=None)[source]
Bases:
objectMohr-Coulomb failure criterion analysis.
Mohr-Coulomb analyzer.
- Parameters:
surface (vtkDataSet) – Surface mesh to analyze with stress data
cohesion (float) – Cohesion in bar
frictionAngle (float) – Friction angle in degrees
logger (Union[Logger, None], optional) – A logger to manage the output messages. Defaults to None, an internal logger is used.
ProfileExtractor
- class geos.processing.tools.ProfileExtractor.ProfileExtractor(logger=None)[source]
Bases:
objectUtility class for extracting profiles along fault surfaces.
- Parameters:
logger (Union[Logger, None], optional) – A logger to manage the output messages. Defaults to None, an internal logger is used.
- extractAdaptiveProfile(centers, values, xStart, yStart, zStart=None, stepSize=20.0, maxSteps=500, cellData=None)[source]
Extract a vertical depth profile with automatic fault detection.
The algorithm adaptively follows a vertical sampling strategy guided by detected fault membership inside the provided cell data. It performs:
Finding the closest starting point to the provided (xStart, yStart, zStart).
- Automatically detecting the target fault using the provided
cellData (e.g., fields like
FaultMaskor any other fault-identifying attribute).
- Automatically detecting the target fault using the provided
Filtering the dataset to keep only cells belonging to that fault.
Splitting the remaining dataset into successive Z-depth slices.
- For each slice, selecting the nearest cell in the XY plane to build the
final vertical profile.
- Parameters:
centers (np.ndarray) – Array of cell centers with shape
(nCells, 3).values (np.ndarray) – Scalar values associated with each cell (shape
(nCells,)).xStart (float or ndarray) – Starting X coordinate.
yStart (float or ndarray) – Starting Y coordinate.
zStart (float or ndarray | None) – Starting Z coordinate. If
None, the method uses the highest point near the provided XY start position.stepSize (float) – Vertical step size used when scanning depth layers. Default is 20.0.
maxSteps (int) – Maximum number of vertical layers to traverse. Default is 500.
cellData (vtkCellData | None) – VTK cell data object containing fields such as
attribute,FaultMask, or other identifiers used to detect and filter the target fault.
- Returns:
depth, profile values, X and Y coordinates of the profile path.
- Return type:
tuple[np.ndarray, np.ndarray, np.ndarray, np.ndarray]
SensitivityAnalyzer
- class geos.processing.tools.SensitivityAnalyzer.SensitivityAnalyzer(outputDir='.', logger=None)[source]
Bases:
objectPerforms sensitivity analysis on Mohr-Coulomb parameters.
Init.
- Parameters:
outputDir (str, optional) – Output directory. Defaults is current directory.
showPlots (bool, optional) – Flag to show the plots. Defaults is True.
logger (Union[Logger, None], optional) – A logger to manage the output messages. Defaults to None, an internal logger is used.
- runAnalysis(surfaceWithStress, time, sensitivityFrictionAngles, sensitivityCohesions, profileStartPoints, profileSearchRadius)[source]
Run sensitivity analysis for multiple friction angles and cohesions.
- Parameters:
surfaceWithStress (vtkDataSet) – Surface to analyze. Should contain stress attribute
time (float) – Time
sensitivityFrictionAngles (list[float]) – List of friction angles to analyze (in degrees)
sensitivityCohesions (list[float]) – List of cohesion to analyze (in bar)
profileStartPoints (list[tuple[float, float]]) – List of start points for profile analysis
profileSearchRadius (float) – Searching radius for determination of profile.
- Returns:
Metrics from input surface.
- Return type:
dict[str, Any]
StressProjector
- class geos.processing.tools.StressProjector.StressProjector(adjacencyMapping, geometricProperties, outputDir='.', logger=None)[source]
Bases:
objectProjects volume stress onto fault surfaces and tracks principal stresses in VTU.
Initialize with pre-computed adjacency mapping and geometric properties.
- Parameters:
adjacencyMapping (dict[int, list[ vtkUnstructuredGrid]]) – Pre-computed dict mapping fault cells to volume cells
geometricProperties (dict[str, Any]) – Pre-computed geometric properties from FaultGeometry: - ‘volumes’: cell volumes - ‘centers’: cell centers - ‘distances’: distances to fault - ‘faultTree’: KDTree for fault
outputDir (str, optional) – Output directory to save plots. Defaults is current directory
logger (Union[Logger, None], optional) – A logger to manage the output messages. Defaults to None, an internal logger is used.
- static computePrincipalStresses(stressTensor)[source]
Compute principal stresses and directions.
Convention: Compression is NEGATIVE - sigma1 = most compressive (most negative) - sigma3 = least compressive (least negative, or most tensile)
- Parameters:
stressTensor (StressTensor) – Stress tensor object
- Returns:
dict with eigenvalues, eigenvectors, meanStress, deviatoricStress
- Return type:
dict[str, npt.NDArray[ np.float64]]
- projectStressToFault(volumeData, volumeInitial, faultSurface, time, timestep=None, weightingScheme=StressProjectorWeightingScheme.ARITHMETIC, computePrincipalStresses=False, frictionAngle=10, cohesion=0)[source]
Project stress and save principal stresses to VTU.
Now uses pre-computed geometric properties for efficiency
- Parameters:
volumeData (vtkUnstructuredGrid) – Volumic mesh.
volumeInitial (vtkUnstructuredGrid) – Pre-simulation mesh
faultSurface (vtkUnstructureGrid) – Fault mesh.
time (float) – Time.
timestep (int | None, optional) – Timestep considered. Defaults is None.
weightingScheme (StressProjectorWeightingScheme, optional) – Weighting scheme for projection. Defaults is “arithmetic”.
computePrincipalStresses (bool, optional) – Flag to compute principal stresses. Defaults is False.
frictionAngle (float, optional) – Friction angle in degrees. Defaults is 10 degrees.
cohesion (float, optional) – Rock cohesion in bar. Defaults is 0 bar.
- Returns:
Fault mesh, volumic mesh, cell contributing mesh.
- Return type:
tuple[vtkUnstructuredGrid, vtkUnstructuredGrid, vtkUnstructuredGrid]
- savePVDCollection(filename='principal_stresses.pvd')[source]
Create PVD file for time series visualization in ParaView.
- Parameters:
filename – Name of PVD file. Defaults is “principal_stresses.pvd”.
- setBiotCoefficientName(biotName)[source]
Set the stress attribute name.
- Parameters:
biotName (str) – Biot coefficient attribute name in the input mesh.