vtk Filters

This package defines vtk filters that allows to process Geos outputs.

geos_posp.filters.GeomechanicsCalculator module

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

GeomechanicsCalculator filter inputs are either vtkPointSet or vtkUnstructuredGrid and returned object is of same type as input.

To use the filter:

import numpy as np
from geos_posp.filters.GeomechanicsCalculator import GeomechanicsCalculator

# filter inputs
logger :Logger
# input object
input :Union[vtkPointSet, vtkUnstructuredGrid]
# grain bulk modulus in Pa (e.g., use a very high value to get a Biot coefficient equal to 1)
grainBulkModulus :float = 1e26
# Reference density to compute specific gravity (e.g. fresh water density) in kg.m^-3
specificDensity :float = 1000.
# rock cohesion in Pa
rockCohesion :float = 1e8
# friction angle in °
frictionAngle :float = 10 * np.pi / 180.

# instanciate the filter
geomechanicsCalculatorFilter :GeomechanicsCalculator = GeomechanicsCalculator()

# set filter attributes
# set logger
geomechanicsCalculatorFilter.SetLogger(logger)
# set input object
geomechanicsCalculatorFilter.SetInputDataObject(input)
# set computeAdvancedOutputsOn or computeAdvancedOutputsOff to compute or
# not advanced outputs...
geomechanicsCalculatorFilter.computeAdvancedOutputsOn()
# set oter parameters
geomechanicsCalculatorFilter.SetGrainBulkModulus(grainBulkModulus)
geomechanicsCalculatorFilter.SetSpecificDensity(specificDensity)
# rock cohesion and friction angle are used for advanced outputs only
geomechanicsCalculatorFilter.SetRockCohesion(rockCohesion)
geomechanicsCalculatorFilter.SetFrictionAngle(frictionAngle)
# do calculations
geomechanicsCalculatorFilter.Update()
# get filter output (same type as input)
output :Union[vtkPointSet, vtkUnstructuredGrid] = geomechanicsCalculatorFilter.GetOutputDataObject(0)
class geos_posp.filters.GeomechanicsCalculator.GeomechanicsCalculator[source]

Bases: VTKPythonAlgorithmBase

VTK Filter to perform Geomechanical output computation.

Input object is either a vtkPointSet or a vtkUntructuredGrid.

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

RequestInformation(request, inInfoVec, outInfoVec)[source]

Inherited from VTKPythonAlgorithmBase::RequestInformation.

Parameters:
  • request (vtkInformation) – request

  • inInfoVec (list[vtkInformationVector]) – input objects

  • outInfoVec (vtkInformationVector) – output objects

Returns:

1 if calculation successfully ended, 0 otherwise.

Return type:

int

SetFrictionAngle(frictionAngle)[source]

Set the friction angle.

Parameters:

frictionAngle (float) – friction angle (rad)

SetGrainBulkModulus(grainBulkModulus)[source]

Set the grain bulk modulus.

Parameters:

grainBulkModulus (float) – grain bulk modulus

SetLogger(logger)[source]

Set the logger.

Parameters:

logger (Logger) – logger

SetRockCohesion(rockCohesion)[source]

Set the rock cohesion.

Parameters:

rockCohesion (float) – rock cohesion

SetSpecificDensity(specificDensity)[source]

Set the specific density.

Parameters:

specificDensity (float) – pecific density

checkMandatoryAttributes(onPoints)[source]

Check that mandatory attributes are present in the mesh.

The mesh must contains either the young Modulus and Poisson’s ratio (m_computeYoungPoisson=False) or the bulk and shear moduli (m_computeYoungPoisson=True)

Parameters:

onPoints (bool) – attributes are on points (True) or on cells (False)

Returns:

True if all needed attributes are present, False otherwise

Return type:

bool

computeAdvancedOutputs()[source]

Compute advanced geomechanical outputs.

Returns:

return True if calculation successfully ended, False otherwise.

Return type:

bool

computeAdvancedOutputsOff()[source]

Deactivate advanced outputs calculation.

computeAdvancedOutputsOn()[source]

Activate advanced outputs calculation.

computeBasicOutputs()[source]

Compute basic geomechanical outputs.

Returns:

return True if calculation successfully ended, False otherwise.

Return type:

bool

computeBiotCoefficient()[source]

Compute Biot coefficient from default and grain bulk modulus.

Returns:

True if calculation successfully ended, False otherwise.

Return type:

bool

computeCompressibilityCoefficient()[source]

Compute compressibility coefficient from simulation outputs.

Compressibility coefficient is computed from Poisson’s ratio, bulk modulus, Biot coefficient and Porosity.

Returns:

True if the attribute is correctly created, False otherwise.

Return type:

bool

computeCriticalPorePressure()[source]

Compute the critical pore pressure and the pressure index.

Returns:

return True if calculation successfully ended, False otherwise.

Return type:

bool

computeCriticalTotalStressRatio()[source]

Compute fracture index and fracture threshold.

Returns:

return True if calculation successfully ended, False otherwise.

Return type:

bool

computeDepthAlongLine()[source]

Compute depth along a line.

Returns:

1D array with depth property

Return type:

npt.NDArray[np.float64]

computeDepthInMesh()[source]

Compute depth of each cell in a mesh.

Returns:

array with depth property

Return type:

npt.NDArray[np.float64]

computeEffectiveStressRatioOed()[source]

Compute the effective stress ratio in oedometric conditions.

Returns:

return True if calculation successfully ended, False otherwise.

Return type:

bool

computeElasticModulus()[source]

Compute elastic moduli.

Young modulus and the poisson ratio computed from shear and bulk moduli if needed.

Returns:

True if elastic moduli are already present or if calculation successfully ended, False otherwise.

Return type:

bool

computeElasticModulusFromBulkShear()[source]

Compute Young modulus Poisson’s ratio from bulk and shear moduli.

Returns:

True if calculation successfully ended, False otherwise

Return type:

bool

computeElasticModulusFromYoungPoisson()[source]

Compute bulk modulus from Young Modulus and Poisson’s ratio.

Returns:

True if bulk modulus was wuccessfully computed, False otherwise

Return type:

bool

computeElasticStrain()[source]

Compute elastic strain from effective stress and elastic modulus.

Returns:

return True if calculation successfully ended, False otherwise.

Return type:

bool

computeLitostaticStress()[source]

Compute lithostatic stress.

Returns:

True if calculation successfully ended, False otherwise.

Return type:

bool

computeRealEffectiveStressRatio()[source]

Compute effective stress ratio.

Returns:

True if calculation successfully ended, False otherwise.

Return type:

bool

computeReservoirStressPathOed()[source]

Compute Reservoir Stress Path in oedometric conditions.

Returns:

return True if calculation successfully ended, False otherwise.

Return type:

bool

computeReservoirStressPathReal()[source]

Compute reservoir stress paths.

Returns:

return True if calculation successfully ended, False otherwise.

Return type:

bool

computeSpecificGravity()[source]

Create Specific gravity attribute.

Specific gravity is computed from rock density attribute and specific density input.

Returns:

True if the attribute is correctly created, False otherwise.

Return type:

bool

computeStressRatioReal(inputAttribute, outputAttribute)[source]

Compute the ratio between horizontal and vertical effective stress.

Returns:

return True if calculation successfully ended, False otherwise.

Return type:

bool

computeTotalStressInitial()[source]

Compute total stress at initial time step.

Returns:

True if calculation successfully ended, False otherwise.

Return type:

bool

computeTotalStresses()[source]

Compute total stress total stress ratio.

Total stress is computed at the initial and current time steps. total stress ratio is computed at current time step only.

Returns:

True if calculation successfully ended, False otherwise.

Return type:

bool

doComputeTotalStress(effectiveStress, pressure, biotCoefficient, totalStressAttributeName)[source]

Compute total stress from effective stress and pressure.

Parameters:
  • effectiveStress (npt.NDArray[np.float64]) – effective stress

  • pressure (npt.NDArray[np.float64] | None) – pressure

  • biotCoefficient (npt.NDArray[np.float64]) – biot coefficient

  • totalStressAttributeName (str) – total stress attribute name

Returns:

return True if calculation successfully ended, False otherwise.

Return type:

bool

getOutputType()[source]

Get output object type.

Returns:

type of output object.

Return type:

str

getPointCoordinates()[source]

Get the coordinates of Points or Cell center.

Returns:

points/cell center coordinates

Return type:

npt.NDArray[np.float64]

getZcoordinates()[source]

Get z coordinates from self.m_output.

Returns:

1D array with depth property

Return type:

npt.NDArray[np.float64]

initFilter()[source]

Check that mandatory attributes are present in the data set.

Determine if attributes are on cells or on Points. Set self.m_ready = True if all data is ok, False otherwise

resetDefaultValues()[source]

Reset filter parameters to the default values.

geos_posp.filters.GeosBlockExtractor module

GeosBlockExtractor module is a vtk filter that allows to extract Volume mesh, Wells and Faults from a vtkMultiBlockDataSet.

Filter input and output types are vtkMultiBlockDataSet.

To use the filter:

from filters.GeosBlockExtractor import GeosBlockExtractor

# filter inputs
logger :Logger
input :vtkMultiBlockDataSet

# instanciate the filter
blockExtractor :GeosBlockExtractor = GeosBlockExtractor()
# set the logger
blockExtractor.SetLogger(logger)
# set input data object
blockExtractor.SetInputDataObject(input)
# ExtractFaultsOn or ExtractFaultsOff to (de)activate the extraction of Fault blocks
blockExtractor.ExtractFaultsOn()
# ExtractWellsOn or ExtractWellsOff to (de)activate the extraction of well blocks
blockExtractor.ExtractWellsOn()
# do calculations
blockExtractor.Update()
# get output object
output :vtkMultiBlockDataSet = blockExtractor.GetOutputDataObject(0)
class geos_posp.filters.GeosBlockExtractor.GeosBlockExtractor[source]

Bases: VTKPythonAlgorithmBase

VTK Filter that perform GEOS block extraction.

The filter returns the volume mesh as the first output, Surface mesh as the second output, and well mesh as the third output.

ExtractFaultsOff()[source]

Deactivate the extraction of Faults.

ExtractFaultsOn()[source]

Activate the extraction of Faults.

ExtractWellsOff()[source]

Deactivate the extraction of Wells.

ExtractWellsOn()[source]

Activate the extraction of Wells.

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

RequestInformation(request, inInfoVec, outInfoVec)[source]

Inherited from VTKPythonAlgorithmBase::RequestInformation.

Parameters:
  • request (vtkInformation) – request

  • inInfoVec (list[vtkInformationVector]) – input objects

  • outInfoVec (vtkInformationVector) – output objects

Returns:

1 if calculation successfully ended, 0 otherwise.

Return type:

int

SetLogger(logger)[source]

Set the logger.

Parameters:

logger (Logger) – logger

UpdateOutputPorts()[source]

Set the number of output ports and update extraction options.

doExtraction()[source]

Apply block extraction.

Returns:

True if block extraction successfully ended, False otherwise.

Return type:

bool

extractRegion(type)[source]

Extract volume mesh from input vtkMultiBlockDataSet.

Returns:

True if volume mesh extraction successfully ended, False otherwise.

Return type:

bool

getOutputFaults()[source]

Get output fault mesh.

Returns:

fault mesh

Return type:

vtkMultiBlockDataSet

getOutputVolume()[source]

Get output volume mesh.

Returns:

volume mesh

Return type:

vtkMultiBlockDataSet

getOutputWells()[source]

Get output well mesh.

Returns:

well mesh

Return type:

vtkMultiBlockDataSet

geos_posp.filters.GeosBlockMerge module

GeosBlockMerge module is a vtk filter that allows to merge Geos ranks, rename stress and porosity attributes, and identify fluids and rock phases.

Filter input and output types are vtkMultiBlockDataSet.

To use the filter:

from filters.GeosBlockMerge import GeosBlockMerge

# filter inputs
logger :Logger
input :vtkMultiBlockDataSet

# instanciate the filter
mergeBlockFilter :GeosBlockMerge = GeosBlockMerge()
# set the logger
mergeBlockFilter.SetLogger(logger)
# set input data object
mergeBlockFilter.SetInputDataObject(input)
# ConvertSurfaceMeshOff or ConvertSurfaceMeshOn to (de)activate the conversion
# of vtkUnstructuredGrid to surface withNormals and Tangents calculation.
mergeBlockFilter.ConvertSurfaceMeshOff()
# do calculations
mergeBlockFilter.Update()
# get output object
output :vtkMultiBlockDataSet = mergeBlockFilter.GetOutputDataObject(0)
class geos_posp.filters.GeosBlockMerge.GeosBlockMerge[source]

Bases: VTKPythonAlgorithmBase

VTK Filter that perform GEOS rank merge.

The filter returns a multiblock mesh composed of elementary blocks.

ConvertSurfaceMeshOff()[source]

Deactivate surface conversion from vtkUnstructredGrid to vtkPolyData.

ConvertSurfaceMeshOn()[source]

Activate surface conversion from vtkUnstructredGrid to vtkPolyData.

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

RequestInformation(request, inInfoVec, outInfoVec)[source]

Inherited from VTKPythonAlgorithmBase::RequestInformation.

Parameters:
  • request (vtkInformation) – request

  • inInfoVec (list[vtkInformationVector]) – input objects

  • outInfoVec (vtkInformationVector) – output objects

Returns:

1 if calculation successfully ended, 0 otherwise.

Return type:

int

SetLogger(logger)[source]

Set the logger.

Parameters:

logger (Logger) – logger

computeNormals(surface)[source]

Compute normals of the given surface.

Parameters:

surface (vtkPolyData) – surface to compute the normals

Returns:

surface with normal attribute

Return type:

vtkPolyData

computeTangents(surface)[source]

Compute tangents of the given surface.

Parameters:

surface (vtkPolyData) – surface to compute the tangents

Returns:

surface with tangent attribute

Return type:

vtkPolyData

convertBlockToSurface(block)[source]

Convert vtkUnstructuredGrid to a surface vtkPolyData.

Warning

work only with triangulated surfaces

Todo

need to convert quadrangular to triangulated surfaces first

Parameters:

block (vtkUnstructuredGrid) – block from which to extract the surface

Returns:

extracted surface

Return type:

vtkPolyData

convertFaultsToSurfaces()[source]

Convert blocks corresponding to faults to surface.

Returns:

True if calculation succesfully ended, False otherwise

Return type:

bool

doMerge()[source]

Apply block merge.

Returns:

True if block merge successfully ended, False otherwise.

Return type:

bool

getPhaseNames(attributeSet)[source]

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

Parameters:

attributeSet (set[str]) – list of attributes where to find phases

Returns:

the list of phase names that appear at least twice.

Return type:

set[str]

getPhases(onCells=True)[source]

Get the dictionnary of phases classified according to PhaseTypeEnum.

Parameters:

onCells (bool, optional) –

Attributes are on Cells (Default) or on Points.

Defaults to True.

Returns:

a dictionnary with phase names as keys and phase type as value.

Return type:

dict[str, PhaseTypeEnum]

mergeChildBlocks(compositeBlock)[source]

Merge all children of the input composite block.

Parameters:

compositeBlock (vtkMultiBlockDataSet) – composite block

Returns:

merged block

Return type:

vtkUnstructuredGrid

mergeRankBlocks()[source]

Merge all elementary node that belong to a same parent node.

Returns:

True if calculation successfully ended, False otherwise

Return type:

bool

renameAttributes(mesh, phaseClassification)[source]

Rename attributes to harmonize throughout the mesh.

Parameters:
  • mesh (vtkMultiBlockDataSet) – input mesh

  • phaseClassification (dict[str, PhaseTypeEnum]) – phase classification detected from attributes

Returns:

output mesh with renamed attributes

Return type:

vtkMultiBlockDataSet

geos_posp.filters.SurfaceGeomechanics module

SurfaceGeomechanics module is a vtk filter that allows to compute Geomechanical properties on surfaces.

Filter input and output types are vtkPolyData.

To use the filter:

from filters.SurfaceGeomechanics import SurfaceGeomechanics

# filter inputs
logger :Logger
input :vtkPolyData
rockCohesion :float
frictionAngle :float # angle in radian

# instanciate the filter
surfaceGeomechanicsFilter :SurfaceGeomechanics = SurfaceGeomechanics()
# set the logger
surfaceGeomechanicsFilter.SetLogger(logger)
# set input data object
surfaceGeomechanicsFilter.SetInputDataObject(input)
# set rock cohesion anf friction angle
surfaceGeomechanicsFilter.SetRockCohesion(rockCohesion)
surfaceGeomechanicsFilter.SetFrictionAngle(frictionAngle)
# do calculations
surfaceGeomechanicsFilter.Update()
# get output object
output :vtkPolyData = surfaceGeomechanicsFilter.GetOutputDataObject(0)
# get created attribute names
newAttributeNames :set[str] = surfaceGeomechanicsFilter.GetNewAttributeNames()
class geos_posp.filters.SurfaceGeomechanics.SurfaceGeomechanics[source]

Bases: VTKPythonAlgorithmBase

Vtk filter to compute geomechanical surfacic attributes.

Input and Output objects are a vtkMultiBlockDataSet containing surface objects with Normals and Tangential attributes.

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

GetNewAttributeNames()[source]

Get the set of new attribute names that were created.

Returns:

set of new attribute names

Return type:

set[str]

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

RequestInformation(request, inInfoVec, outInfoVec)[source]

Inherited from VTKPythonAlgorithmBase::RequestInformation.

Parameters:
  • request (vtkInformation) – request

  • inInfoVec (list[vtkInformationVector]) – input objects

  • outInfoVec (vtkInformationVector) – output objects

Returns:

1 if calculation successfully ended, 0 otherwise.

Return type:

int

SetFrictionAngle(frictionAngle)[source]

Set friction angle value.

Parameters:

frictionAngle (float) – friction angle (rad)

SetLogger(logger)[source]

Set the logger.

Parameters:

logger (Logger) – logger

SetRockCohesion(rockCohesion)[source]

Set rock cohesion value.

Parameters:

rockCohesion (float) – rock cohesion (Pa)

computeChangeOfBasisMatrix(localBasis, fromLocalToYXZ)[source]

Compute the change of basis matrix according to local coordinates.

Parameters:
  • localBasis (npt.NDArray[np.float64]) – local coordinate vectors.

  • fromLocalToYXZ (bool) – if True, change of basis matrix is from local to XYZ bases, otherwise it is from XYZ to local bases.

Returns:

change of basis matrix.

Return type:

npt.NDArray[np.float64]

computeNewCoordinates(attrArray, normalTangentVectors, fromLocalToYXZ)[source]

Compute the coordinates of a vectorial attribute.

Parameters:
  • attrArray (npt.NDArray[np.float64]) – vector of attribute values

  • normalTangentVectors (npt.NDArray[np.float64]) – 3xNx3 local vector coordinates

  • fromLocalToYXZ (bool) – if True, conversion is done from local to XYZ basis, otherwise conversion is done from XZY to Local basis.

Returns:

Vector of new coordinates of the attribute.

Return type:

npt.NDArray[np.float64]

computeNewCoordinatesVector3(vector, changeOfBasisMatrix)[source]

Compute attribute new coordinates of vector of size 3.

Parameters:
  • vector (npt.NDArray[np.float64]) – input coordinates.

  • changeOfBasisMatrix (npt.NDArray[np.float64]) – change of basis matrix

Returns:

new coordinates

Return type:

npt.NDArray[np.float64]

computeNewCoordinatesVector6(vector, changeOfBasisMatrix)[source]

Compute attribute new coordinates of vector of size > 3.

Parameters:
  • vector (npt.NDArray[np.float64]) – input coordinates.

  • changeOfBasisMatrix (npt.NDArray[np.float64]) – change of basis matrix

Returns:

new coordinates

Return type:

npt.NDArray[np.float64]

computeShearCapacityUtilization()[source]

Compute the shear capacity utilization on surface.

Returns:

True if calculation successfully ended, False otherwise.

Return type:

bool

convertLocalToXYZBasisAttributes()[source]

Convert vectorial property coordinates from Local to canonic basis.

Returns:

True if calculation successfully ended, False otherwise

Return type:

bool

convertXYZToLocalBasisAttributes()[source]

Convert vectorial property coordinates from canonic to local basis.

Returns:

True if calculation successfully ended, False otherwise

Return type:

bool

filterAttributesToConvert(attributesToFilter0)[source]

Filter the set of attribute names if they are vectorial and present.

Parameters:

attributesToFilter0 (set[str]) – set of attribute names to filter.

Returns:

Set of the attribute names.

Return type:

set[str]

getAttributesToConvertFromLocalToXYZ()[source]

Get the list of attributes to convert from local to XYZ basis.

Returns:

Set of the attribute names.

Return type:

set[str]

getNormalTangentsVectors()[source]

Compute the change of basis matrix from Local to XYZ bases.

Returns:

Nx3 matrix of local vector coordinates.

Return type:

npt.NDArray[np.float64]