Post-processing filters

This module contains GEOS post-processing tools.

Geomechanics post-processing tools

GEOS computes many outputs including flow and geomechanic properties if coupled flow geomechanics simulations were run. Users however need additional metrics to asses geomechanical stability as part of operational studies. Two filters have been developped to compute these additional metrics: Geos Geomechanics Calculator` and Geos Surfacic Geomechanics. In addition, a third filter allows to plot Mohr’s circles on selected cells and time steps.

Note

Several processing filters require the definition of physical parameters. The following list is non-exhaustive. See PhysicalConstants module.

Default values: - grainBulkModulus = 38e9 Pa ( quartz value ) - specificDensity = 1000.0 kg/m³ ( water value ) - rock cohesion: 0.0 Pa $( fractured case ) - friction angle: 10°

SurfaceGeomechanics

SurfaceGeomechanics is a VTK filter that allows:
  • Conversion of a set of attributes from local basis to XYZ basis

  • Computation of the shear capacity utilization (SCU)

Warning

  • The computation of the SCU requires the presence of a ‘traction’ attribute in the input mesh.

  • Conversion from local to XYZ basis is currently only handled for cell attributes.

Note

Default values for physical constants used in this filter:
  • rock cohesion: 0.0 Pa ( fractured case)

  • friction angle: 10°

Filter input and output types are vtkPolyData.

To use the filter:

from geos.processing.post_processing.SurfaceGeomechanics import SurfaceGeomechanics
from geos.utils.Errors import VTKError

# filter inputs
inputMesh: vtkPolyData

# Optional inputs
rockCohesion: float   # Pa
frictionAngle: float  # angle in radian
speHandler: bool      # defaults to false

# Instantiate the filter
sg: SurfaceGeomechanics = SurfaceGeomechanics( inputMesh, speHandler )

# Set the handler of yours (only if speHandler is True).
yourHandler: logging.Handler
sg.SetLoggerHandler( yourHandler )

# [optional] Set rock cohesion and friction angle
sg.SetRockCohesion( rockCohesion )
sg.SetFrictionAngle( frictionAngle )

# Do calculations
try:
    sg.applyFilter()
except ( ValueError, VTKError, AttributeError, AssertionError ) as e:
    sg.logger.error( f"The filter { sg.logger.name } failed due to: { e }" )
except Exception as e:
    mess: str = f"The filter { sg.logger.name } failed due to: { e }"
    sg.logger.critical( mess, exc_info=True )

# Get output object
output: vtkPolyData = sg.GetOutputMesh()

# Get created attribute names
newAttributeNames: set[str] = sg.GetNewAttributeNames()

Note

By default, conversion of attributes from local to XYZ basis is performed for the following list: { ‘displacementJump’ }. This list can be modified in different ways:

  • Addition of one or several additional attributes to the set by using the filter function AddAttributesToConvert.

  • Replace the list completely with the function SetAttributesToConvert.

Note that the dimension of the attributes to convert must be equal or greater than 3.

class geos.processing.post_processing.SurfaceGeomechanics.SurfaceGeomechanics(surfacicMesh, speHandler=False)[source]

Bases: object

Vtk filter to compute geomechanical surfacic attributes.

Input and Output objects are a vtkPolydata with surfaces objects with Normals and Tangential attributes.

Parameters:
  • surfacicMesh (vtkPolyData) – The input surfacic mesh.

  • speHandler (bool, optional) – True to use a specific handler, False to use the internal handler. Defaults to False.

AddAttributesToConvert(attributeName)[source]

Add an attribute to the set of attributes to convert.

Parameters:

attributeName (Union[list[str],set[str]]) – List of the attribute array names.

ConvertAttributesOff()[source]

Deactivate the conversion of attributes from local to XYZ basis.

ConvertAttributesOn()[source]

Activate the conversion of attributes from local to XYZ basis.

GetAttributesToConvert()[source]

Get the set of attributes that will be converted from local to XYZ basis.

Returns:

The set of attributes that should be converted.

Return type:

set[ str ]

GetConvertAttributes()[source]

If convertAttributesOn is True, the set of attributes will be converted from local to XYZ during filter application. Default is True.

Returns:

Current state of the attributes conversion request.

Return type:

bool

GetNewAttributeNames()[source]

Get the set of new attribute names that were created.

Returns:

The set of new attribute names.

Return type:

set[str]

GetOutputMesh()[source]

Get the output mesh with computed attributes.

Returns:

The surfacic output mesh.

Return type:

vtkPolyData

SetAttributesToConvert(attributesToConvert)[source]

Set the list of attributes that will be converted from local to XYZ basis.

Parameters:

attributesToConvert (set[str]) – The set of attributes names that will be converted from local to XYZ basis

SetFrictionAngle(frictionAngle)[source]

Set friction angle value in radians. Defaults to 10 / 180 * pi rad.

Parameters:

frictionAngle (float) – The friction angle in radians.

SetLoggerHandler(handler)[source]

Set a specific handler for the filter logger.

In this filter 4 log levels are use, .info, .error, .warning and .critical, be sure to have at least the same 4 levels.

Parameters:

handler (logging.Handler) – The handler to add.

SetRockCohesion(rockCohesion)[source]

Set rock cohesion value. Defaults to 0.0 Pa.

Parameters:

rockCohesion (float) – The rock cohesion in Pascal.

SetSurfaceName(name)[source]

Set a name for the input surface. For logging purpose only.

Parameters:

name (str) – The identifier for the surface.

applyFilter()[source]

Compute Geomechanical properties on input surface.

Raises:
  • ValueError – Errors during the creation of an attribute.

  • VTKError – Error raises during the call of VTK function.

  • AttributeError – Attributes must be on cell.

  • AssertionError – Something went wrong during the shearCapacityUtilization computation.

computeShearCapacityUtilization()[source]

Compute the shear capacity utilization (SCU) on surface.

Raises:
  • ValueError – Something went wrong during the creation of an attribute.

  • AssertionError – Something went wrong during the shearCapacityUtilization computation.

convertAttributesFromLocalToXYZBasis()[source]

Convert attributes from local to XYZ basis.

Raises:
  • ValueError – Something went wrong during the creation of an attribute.

  • AttributeError – Attributes must be on cell.

GeomechanicsCalculator

GeomechanicsCalculator module is a vtk filter that allows to compute additional geomechanical properties from existing ones on a vtkUnstructuredGrid.

To compute the geomechanics outputs, the mesh must have the following properties:
  • The Young modulus and the Poisson’s ratio named “youngModulus” and “poissonRatio” or bulk and shear moduli named “bulkModulus” and “shearModulus”

  • The initial Young modulus and Poisson’s ratio named “youngModulusInitial” and “poissonRatioInitial” or the initial bulk modulus named “bulkModulusInitial”

  • The porosity named “porosity”

  • The initial porosity named “porosityInitial”

  • The delta of pressure named “deltaPressure”

  • The density named “density”

  • The effective stress named “stressEffective”

  • The initial effective stress named “stressEffectiveInitial”

  • The pressure named “pressure”

The basic geomechanics properties computed on the mesh are:
  • The elastic moduli not present on the mesh

  • Biot coefficient

  • Compressibility, oedometric compressibility and real compressibility coefficient

  • Specific gravity

  • Real effective stress ratio

  • Total initial stress, total current stress and total stress ratio

  • Elastic stain

  • Real reservoir stress path and reservoir stress path in oedometric condition

The advanced geomechanics properties computed on the mesh are:
  • Fracture index and threshold

  • Critical pore pressure and pressure index

The input and output meshes are vtkUnstructuredGrid

Note

  • The default physical constants used by the filter are:
    • grainBulkModulus = 38e9 Pa ( quartz value )

    • specificDensity = 1000.0 kg/m³ ( water value )

    • rockCohesion = 0.0 Pa ( fractured case )

    • frictionAngle = 10.0°

To use the filter:

import numpy as np
from geos.processing.post_processing.GeomechanicsCalculator import GeomechanicsCalculator

# Define filter inputs
mesh: vtkUnstructuredGrid
computeAdvancedProperties: bool # optional, defaults to False
speHandler: bool # optional, defaults to False
loggerName: str # Defaults to "Geomechanics Calculator"

# Instantiate the filter
geomechanicsCalculatorFilter: GeomechanicsCalculator = GeomechanicsCalculator( mesh, computeAdvancedProperties, speHandler )

# Use your own handler (if speHandler is True)
yourHandler: logging.Handler
geomechanicsCalculatorFilter.setLoggerHandler( yourHandler )

# Change the physical constants if needed
## For the basic properties
grainBulkModulus: float
specificDensity: float
geomechanicsCalculatorFilter.physicalConstants.grainBulkModulus = grainBulkModulus
geomechanicsCalculatorFilter.physicalConstants.specificDensity = specificDensity

## For the advanced properties
rockCohesion: float
frictionAngle: float
geomechanicsCalculatorFilter.physicalConstants.rockCohesion = rockCohesion
geomechanicsCalculatorFilter.physicalConstants.frictionAngle = frictionAngle

# Do calculations
try:
    geomechanicsCalculatorFilter.applyFilter()
except ( ValueError, AttributeError ) as e:
    geomechanicsCalculatorFilter.logger.error( f"The filter { geomechanicsCalculatorFilter.logger.name } failed due to: { e }" )
except Exception as e:
    mess: str = f"The filter { geomechanicsCalculatorFilter.logger.name } failed due to: { e }"
    geomechanicsCalculatorFilter.logger.critical( mess, exc_info=True )

# Get the mesh with the geomechanics properties computed as attribute
output: vtkUnstructuredGrid
output = geomechanicsCalculatorFilter.getOutput()
class geos.processing.post_processing.GeomechanicsCalculator.GeomechanicsCalculator(mesh, computeAdvancedProperties=False, loggerName='Geomechanics Calculator', speHandler=False)[source]

Bases: object

VTK Filter to perform geomechanics properties computation.

Parameters:
  • mesh (vtkUnstructuredGrid) – Input mesh.

  • computeAdvancedProperties (bool, optional) – True to compute advanced geomechanics properties, False otherwise. Defaults to False.

  • loggerName (str, optional) – Name of the filter logger. Defaults to “Geomechanics Calculator”.

  • speHandler (bool, optional) – True to use a specific handler, False to use the internal handler. Defaults to False.

class AdvancedProperties(_criticalTotalStressRatio=None, _stressRatioThreshold=None, _criticalPorePressure=None, _criticalPorePressureIndex=None)[source]

Bases: object

The dataclass with the value of the advanced geomechanics properties.

property criticalPorePressure

Get the critical pore pressure value.

property criticalPorePressureIndex

Get the critical pore pressure index value.

property criticalTotalStressRatio

Get the critical total stress ratio value.

getAdvancedPropertyValue(name)[source]

Get the value of the advanced property wanted.

Parameters:

name (str) – The name of the advanced property.

Returns:

the value of the advanced property.

Return type:

npt.NDArray[ np.float64 ]

property stressRatioThreshold

Get the stress ratio threshold value.

class BasicProperties(_biotCoefficient=None, _compressibility=None, _compressibilityOed=None, _compressibilityReal=None, _specificGravity=None, _effectiveStressRatioReal=None, _totalStress=None, _totalStressT0=None, _totalStressRatioReal=None, _elasticStrain=None, _deltaTotalStress=None, _rspReal=None, _rspOed=None, _effectiveStressRatioOed=None)[source]

Bases: object

The dataclass with the value of the basic geomechanics properties.

property biotCoefficient

Get the biot coefficient value.

property compressibility

Get the compressibility value.

property compressibilityOed

Get the compressibility in oedometric condition value.

property compressibilityReal

Get the real compressibility value.

property deltaTotalStress

Get the total delta stress value.

property effectiveStressRatioOed

Get the effective stress ratio oedometric value.

property effectiveStressRatioReal

Get the real effective stress ratio value.

property elasticStrain

Get the elastic strain value.

getBasicPropertyValue(name)[source]

Get the value of the basic property wanted.

Parameters:

name (str) – The name of the basic property.

Returns:

the value of the basic property.

Return type:

npt.NDArray[ np.float64 ]

property rspOed

Get the reservoir stress path in oedometric condition value.

property rspReal

Get the real reservoir stress path value.

property specificGravity

Get the specific gravity value.

property totalStress

Get the total stress value.

property totalStressRatioReal

Get the total real stress ratio value.

property totalStressT0

Get the initial total stress value.

class ElasticModuli(_bulkModulus=None, _shearModulus=None, _youngModulus=None, _poissonRatio=None, _bulkModulusT0=None, _youngModulusT0=None, _poissonRatioT0=None)[source]

Bases: object

The dataclass with the value of the elastic moduli.

property bulkModulus

Get the bulk modulus value.

property bulkModulusT0

Get the bulk modulus at the initial time value.

getElasticModulusValue(name)[source]

Get the wanted elastic modulus value.

Parameters:

name (str) – The name of the wanted elastic modulus.

Returns:

The value of the elastic modulus.

Return type:

npt.NDArray[np.float64]

property poissonRatio

Get the poisson ratio value.

property poissonRatioT0

Get the poisson ration at the initial time value.

setElasticModulusValue(name, value)[source]

Set the elastic modulus value wanted.

Parameters:
  • name (str) – The name of the elastic modulus.

  • value (npt.NDArray[np.float64]) – The value to set.

property shearModulus

Get the shear modulus value.

property youngModulus

Get the young modulus value.

property youngModulusT0

Get the young modulus at the initial time value.

class MandatoryProperties(_porosity=None, _porosityInitial=None, _pressure=None, _deltaPressure=None, _density=None, _effectiveStress=None, _effectiveStressT0=None)[source]

Bases: object

The dataclass with the value of mandatory properties to compute the other geomechanics properties.

property deltaPressure

Get the delta pressure value.

property density

Get the density value.

property effectiveStress

Get the effective stress value.

property effectiveStressT0

Get the initial effective stress value.

property porosity

Get the porosity value.

property porosityInitial

Get the initial porosity value.

property pressure

Get the pressure value.

setMandatoryPropertyValue(name, value)[source]

Set the value of a mandatory property.

Parameters:
  • name (str) – The name of the property.

  • value (npt.NDArray[np.float64]) – The value to set.

class PhysicalConstants(_grainBulkModulus=38000000000.0, _specificDensity=1000.0, _rockCohesion=0.0, _frictionAngle=0.17453292519943295)[source]

Bases: object

The dataclass with the value of the physical constant used to compute geomechanics properties.

property frictionAngle

Get the friction angle value.

property grainBulkModulus

Get the grain bulk modulus value.

property rockCohesion

Get the rock cohesion value.

property specificDensity

Get the specific density value.

applyFilter()[source]

Compute the geomechanics properties and create attributes on the mesh.

Raises:
  • AttributeError – A mandatory attribute is missing.

  • ValueError – Something went wrong during the creation of an attribute.

getOutput()[source]

Get the mesh with the geomechanics properties computed as attributes.

Returns:

The mesh with the geomechanics properties computed as attributes.

Return type:

vtkUnstructuredGrid

getOutputType()[source]

Get output object type.

Returns:

Type of output object.

Return type:

str

physicalConstants
setLoggerHandler(handler)[source]

Set a specific handler for the filter logger.

In this filter 4 log levels are use, .info, .error, .warning and .critical, be sure to have at least the same 4 levels.

Parameters:

handler (logging.Handler) – The handler to add.

GeosBlockExtractor

GeosBlockExtractor is a vtk filter that allows to extract the domain (volume, fault and well) from a GEOS output multiBlockDataset mesh.

Important

The input mesh must be an output of a GEOS simulation or contain at least three blocks labeled with the same domain names:
  • “CellElementRegion” for volume domain

  • “SurfaceElementRegion” for fault domain

  • “WellElementRegion” for well domain

See more in GEOS documentation about Element region.

Note

Volume domain is automatically extracted, by defaults Fault and Well domains are empty multiBlockDataSet.

To use the filter:

from geos.processing.post_processing.GeosBlockExtractor import GeosBlockExtractor

# Filter inputs.
geosMesh: vtkMultiBlockDataSet
# Optional inputs.
extractFault: bool # Defaults to False
extractWell: bool # Defaults to False
speHandler: bool # Defaults to False
loggerName: str # Defaults to "Geos Block Extractor"

# Instantiate the filter
geosBlockExtractor: GeosBlockExtractor = GeosBlockExtractor( geosMesh, extractFault, extractWell, speHandler )

# Set the handler of yours (only if speHandler is True).
yourHandler: logging.Handler
geosBlockExtractor.setLoggerHandler( yourHandler )

# Do calculations
try:
    geosBlockExtractor.applyFilter()
except ( ValueError, TypeError ) as e:
    geosBlockExtractor.logger.error( f"The filter { geosBlockExtractor.logger.name } failed due to: { e }." )
except Exception as e:
    mess: str = f"The filter { geosBlockExtractor.logger.name } failed du to: { e }"
    geosBlockExtractor.logger.critical( mess, exc_info=True )

# Get the multiBlockDataSet with blocks of the extracted domain.
geosDomainExtracted: vtkMultiBlockDataSet
geosDomainExtracted = geosBlockExtractor.extractedGeosDomain.volume # For volume domain
geosDomainExtracted = geosBlockExtractor.extractedGeosDomain.fault # For fault domain
geosDomainExtracted = geosBlockExtractor.extractedGeosDomain.well # For well domain
class geos.processing.post_processing.GeosBlockExtractor.GeosBlockExtractor(geosMesh, extractFault=False, extractWell=False, speHandler=False, loggerName='Geos Block Extractor')[source]

Bases: object

Blocks from the ElementRegions from a GEOS output multiBlockDataset mesh.

Parameters:
  • geosMesh (vtkMultiBlockDataSet) – The mesh from Geos.

  • extractFault (bool, optional) – True if SurfaceElementRegion needs to be extracted, False otherwise. Defaults to False.

  • extractWell (bool, optional) – True if WellElementRegion needs to be extracted, False otherwise. Defaults to False.

  • speHandler (bool, optional) – True to use a specific handler, False to use the internal handler. Defaults to False.

  • loggerName (str, optional) – Name of the filter logger. Defaults to “Geos Block Extractor”.

class ExtractedGeosDomain(_volume=<vtkmodules.vtkCommonDataModel.vtkMultiBlockDataSet(0x5a78597f68d0)>, _fault=<vtkmodules.vtkCommonDataModel.vtkMultiBlockDataSet(0x5a78594a9ef0)>, _well=<vtkmodules.vtkCommonDataModel.vtkMultiBlockDataSet(0x5a7859330a50)>)[source]

Bases: object

The dataclass with the three GEOS domain mesh.

property fault

Get the mesh with the blocks of the GEOS SurfaceElementRegion.

setExtractedDomain(geosDomainName, multiBlockDataSet)[source]

Set the mesh to the correct domain.

Parameters:
  • geosDomainName (GeosDomainNameEnum) – Name of the GEOS domain.

  • multiBlockDataSet (vtkMultiBlockDataSet) – The mesh to set.

Raises:

ValueError – The mesh is not a GEOS domain.

property volume

Get the mesh with the blocks of the GEOS CellElementRegion.

property well

Get the mesh with the blocks of the GEOS WellElementRegion.

applyFilter()[source]

Extract the volume, the fault or the well domain of the mesh from GEOS.

Raises:
  • ValueError – The mesh extracted is not a GEOS domain.

  • TypeError – The mesh extracted has the wrong dimension.

extractedGeosDomain
setLoggerHandler(handler)[source]

Set a specific handler for the filter logger.

In this filter 4 log levels are use, .info, .error, .warning and .critical, be sure to have at least the same 4 levels.

Parameters:

handler (logging.Handler) – The handler to add.

class geos.processing.post_processing.GeosBlockExtractor.GeosExtractDomainBlock[source]

Bases: vtkExtractBlock

Extract blocks from a GEOS output multiBlockDataset mesh.

AddGeosDomainName(geosDomainName)[source]

Add the index of the GEOS domain to extract from its name.

Parameters:

geosDomainName (GeosDomainNameEnum) – Name of the GEOS domain to extract.

GeosBlockMerge

GeosBlockMerge is a vtk filter that acts on each region of a GEOS output domain (volume, fault, wells):
  • Ranks are merged.

  • “Fluids” and “Rock” phases are identified.

  • “Rock” attributes are renamed depending on the phase they refer to for more clarity.

  • Volume meshes are converted to surface if needed.

Filter input and output types are vtkMultiBlockDataSet.

Important

This filter cannot be used directly on GEOS output. The domain needs to be extracted with the help of GeosBlockExtractor filter. Please refer to the documentation for more information.

To use the filter:

from geos.processing.post_processing.GeosBlockMerge import GeosBlockMerge
from geos.utils.Errors import VTKError

# Filter inputs.
inputMesh: vtkMultiBlockDataSet
# Optional inputs.
convertFaultToSurface: bool # Defaults to False
speHandler: bool # Defaults to False
loggerName: str # Defaults to "GEOS Block Merge"

# Instantiate the filter
mergeBlockFilter: GeosBlockMerge = GeosBlockMerge( inputMesh, convertFaultToSurface, speHandler, loggerName )

# Set the handler of yours (only if speHandler is True).
yourHandler: logging.Handler
mergeBlockFilter.setLoggerHandler( yourHandler )

# Do calculations
try:
    mergeBlockFilter.applyFilter()
except ( ValueError, VTKError ) as e:
    mergeBlockFilter.logger.error( f"The filter { mergeBlockFilter.logger.name } failed due to: { e }" )
except Exception as e:
    mess: str = f"The filter { mergeBlockFilter.logger.name } failed due to: { e }"
    mergeBlockFilter.logger.critical( mess, exc_info=True )

# Get the multiBlockDataSet with one dataSet per region
outputMesh: vtkMultiBlockDataSet = mergeBlockFilter.getOutput()
class geos.processing.post_processing.GeosBlockMerge.GeosBlockMerge(inputMesh, convertFaultToSurface=False, speHandler=False, loggerName='GEOS Block Merge')[source]

Bases: object

VTK Filter that merges ranks of GEOS output mesh.

For each composite block of the input mesh:
  • Ranks are merged.

  • “Rock” attributes are renamed.

  • Volume meshes are converted to surface if requested.

Parameters:
  • inputMesh (vtkMultiBlockDataSet) – The mesh with the blocks to merge.

  • convertFaultToSurface (bool, optional) – If True, merged blocks are converted to surface (vtp), False otherwise. Defaults to False.

  • speHandler (bool, optional) – True to use a specific handler, False to use the internal handler. Defaults to False.

  • loggerName (str, optional) – Name of the filter logger. Defaults to “GEOS Block Merge”.

applyFilter()[source]

Apply the filter on the mesh.

Raises:
  • ValueError – Something went wrong during the creation of an attribute.

  • VTKError – Error raises during the call of VTK function.

computePhaseNames()[source]

Get the names of the phases in the mesh from Cell attributes.

getOutput()[source]

Get the mesh with the composite blocks merged.

renameAttributes(mesh)[source]

Rename attributes to harmonize GEOS output, see more geos.utils.OutputsConstants.py.

Parameters:

mesh (vtkUnstructuredGrid) – The mesh with the attribute to rename.

setLoggerHandler(handler)[source]

Set a specific handler for the filter logger.

In this filter 4 log levels are use, .info, .error, .warning and .critical, be sure to have at least the same 4 levels.

Parameters:

handler (logging.Handler) – The handler to add.