Processing functions

This package define functions to process data.

geos_posp.processing.geomechanicsCalculatorFunctions module

Functions to compute additional geomechanical properties.

geos_posp.processing.geomechanicsCalculatorFunctions.biotCoefficient(Kgrain, bulkModulus)[source]

Compute Biot coefficient.

\[b = 1-\frac{K}{K_{grain}}\]
Parameters:
  • Kgrain (float) – grain bulk modulus (\(K_{grain}\) - Pa)

  • bulkModulus (npt.NDArray[np.float64]) – default bulk modulus (K - Pa)

Returns:

biot coefficient (b)

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.bulkModulus(youngModulus, poissonRatio)[source]

Compute bulk Modulus from young modulus and poisson ratio.

\[K = \frac{E}{3(1-2\nu)}\]
Parameters:
  • youngModulus (npt.NDArray[np.float64]) – Young modulus (E - Pa)

  • poissonRatio (npt.NDArray[np.float64]) – Poisson’s ratio (\(\nu\))

Returns:

Bulk modulus (K - Pa)

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.compressibility(poissonRatio, bulkModulus, biotCoefficient, porosity)[source]

Compute compressibility from elastic moduli, biot coefficient and porosity.

Compressibility formula is:

\[C = \frac{1-2\nu}{\phi K}\left(\frac{b²(1+\nu)}{1-\nu} + 3(b-\phi)(1-b)\right)\]
Parameters:
  • poissonRatio (npt.NDArray[np.float64]) – Poisson’s ratio (\(\nu\)).

  • bulkModulus (npt.NDArray[np.float64]) – Bulk Modulus (K - Pa)

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

  • porosity (npt.NDArray[np.float64]) – Rock porosity (\(\phi\)).

Returns:

Compressibility array (C - Pa^-1).

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.compressibilityOed(shearModulus, bulkModulus, porosity)[source]

Compute compressibility from elastic moduli and porosity.

Compressibility formula is:

\[C = \frac{1}{\phi}.\frac{3}{3K+4G}\]
Parameters:
  • shearModulus (npt.NDArray[np.float64]) – shear modulus (G - Pa).

  • bulkModulus (npt.NDArray[np.float64]) – Bulk Modulus (K - Pa).

  • porosity (npt.NDArray[np.float64]) – Rock porosity (\(\phi\)).

Returns:

Oedometric Compressibility (C - Pa^-1).

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.compressibilityReal(deltaPressure, porosity, porosityInitial)[source]

Compute compressibility from elastic moduli and porosity.

Compressibility formula is:

\[C = \frac{\phi-\phi_0}{\Delta P.\phi_0}\]
Parameters:
  • deltaPressure (npt.NDArray[np.float64]) – Pressure deviation (\(\Delta P\) - Pa).

  • porosity (npt.NDArray[np.float64]) – Rock porosity (\(\phi\)).

  • porosityInitial (npt.NDArray[np.float64]) – initial porosity (\(\phi_0\)).

Returns:

Real compressibility (C - Pa^-1).

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.computeNormalShearStress(stressTensor, directionVector)[source]

Compute normal and shear stress according to stress tensor and direction.

Parameters:
  • stressTensor (npt.NDArray[np.float64]) – 3x3 stress tensor

  • directionVector (npt.NDArray[np.float64]) – direction vector

Returns:

normal and shear stresses.

Return type:

tuple[float, float]

geos_posp.processing.geomechanicsCalculatorFunctions.computeStressPrincipalComponents(stressTensor)[source]

Compute stress principal components.

Parameters:

stressTensor (npt.NDArray[np.float64]) – stress tensor.

Returns:

Principal components sorted in ascending

order.

Return type:

tuple[float, float, float]

geos_posp.processing.geomechanicsCalculatorFunctions.computeStressPrincipalComponentsFromStressVector(stressVector)[source]

Compute stress principal components from stress vector.

Parameters:

stressVector (npt.NDArray[np.float64]) – stress vector.

Returns:

Principal components sorted in ascending

order.

Return type:

tuple[float, float, float]

geos_posp.processing.geomechanicsCalculatorFunctions.criticalPorePressure(stressVector, rockCohesion, frictionAngle=0.0)[source]

Compute the critical pore pressure.

Fracturing can occur in areas where Critical pore pressure is greater than the pressure. (see Khan, S., Khulief, Y., Juanes, R., Bashmal, S., Usman, M., & Al-Shuhail, A. (2024). Geomechanical Modeling of CO2 Sequestration: A Review Focused on CO2 Injection and Monitoring. Journal of Environmental Chemical Engineering, 112847. https://doi.org/10.1016/j.jece.2024.112847)

\[P_{Cr}=\frac{c.\cos(\alpha)}{1-\sin(\alpha)} + \frac{3\sigma_3-\sigma_1}{2}\]
Parameters:
  • stressVector (npt.NDArray[np.float64]) – stress vector (\(\sigma\) - Pa).

  • rockCohesion (float) – rock cohesion (c - Pa)

  • frictionAngle (float, optional) –

    friction angle (\(\alpha\) - rad).

    Defaults to 0 rad.

Returns:

critical pore pressure (\(P_{Cr}\) - Pa)

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.criticalPorePressureThreshold(pressure, criticalPorePressure)[source]

Compute the critical pore pressure threshold.

Defined as the ratio between pressure and critical pore pressure. Fracturing can occur in areas where critical pore pressure threshold is >1.

\[P_{Th}=\frac{P}{P_{Cr}}\]
Parameters:
  • pressure (npt.NDArray[np.float64]) – pressure (P - Pa)

  • criticalPorePressure (npt.NDArray[np.float64]) – Critical pore pressure (\(P_{Cr}\) - Pa)

Returns:

Critical pore pressure threshold (\(P_{Th}\))

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.criticalTotalStressRatio(pressure, verticalStress)[source]

Compute critical total stress ratio.

Corresponds to the fracture index from Lemgruber-Traby et al (2024). Fracturing can occur in areas where K > Total stress ratio. (see Lemgruber-Traby, A., Cacas, M. C., Bonte, D., Rudkiewicz, J. L., Gout, C., & Cornu, T. (2024). Basin modelling workflow applied to the screening of deep aquifers for potential CO2 storage. Geoenergy, geoenergy2024-010. https://doi.org/10.1144/geoenergy2024-010)

\[\sigma_{Ĉr}=\frac{P}{\sigma_v}\]
Parameters:
  • pressure (npt.NDArray[np.float64]) – Pressure (P - Pa)

  • verticalStress (npt.NDArray[np.float64]) – Vertical total stress (\(\sigma_v\) - Pa) using Geos convention

Returns:

critical total stress ratio

(\(\sigma_{Cr}\))

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.deviatoricStressPathOed(poissonRatio)[source]

Compute the Deviatoric Stress Path parameter in oedometric conditions.

This parameter corresponds to the ratio between horizontal and vertical stresses in oedometric conditions.

\[DSP=\frac{\nu}{1-\nu}\]
Parameters:

poissonRatio (npt.NDArray[np.float64]) – Poisson’s ratio (\(\nu\))

Returns:

Deviatoric Stress Path parameter in

oedometric conditions (DSP)

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.elasticStrainFromBulkShear(deltaEffectiveStress, bulkModulus, shearModulus)[source]

Compute elastic strain from Bulk and Shear moduli.

See documentation on https://dnicolasespinoza.github.io/node5.html.

\[ \begin{align}\begin{aligned}\epsilon=\Delta\sigma_{eff}.C^{-1}\\\begin{split}C=\begin{pmatrix} K+\frac{4}{3}G & K-\frac{2}{3}G & K-\frac{2}{3}G & 0 & 0 & 0\\ K-\frac{2}{3}G & K+\frac{4}{3}G & K-\frac{2}{3}G & 0 & 0 & 0\\ K-\frac{2}{3}G & K-\frac{2}{3}G & K+\frac{4}{3}G & 0 & 0 & 0\\ 0 & 0 & 0 & \nu & 0 & 0\\ 0 & 0 & 0 & 0 & \nu & 0\\ 0 & 0 & 0 & 0 & 0 & \nu\\ \end{pmatrix}\end{split}\end{aligned}\end{align} \]

where C is stiffness tensor.

Parameters:
  • deltaEffectiveStress (npt.NDArray[np.float64]) – effective stress variation (\(\Delta\sigma_{eff}\) - Pa) [S11, S22, S33, S23, S13, S12]

  • bulkModulus (npt.NDArray[np.float64]) – Bulk modulus (K - Pa)

  • shearModulus (npt.NDArray[np.float64]) – Shear modulus (G - Pa)

Returns:

elastic strain (\(\epsilon\))

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.elasticStrainFromYoungPoisson(deltaEffectiveStress, youngModulus, poissonRatio)[source]

Compute elastic strain from Young modulus and Poisson ratio.

See documentation on https://dnicolasespinoza.github.io/node5.html.

\[ \begin{align}\begin{aligned}\epsilon=\Delta\sigma_{eff}.C^{-1}\\\begin{split}C=\begin{pmatrix} \lambda+2G & \lambda & \lambda & 0 & 0 & 0\\ \lambda & \lambda+2G & \lambda & 0 & 0 & 0\\ \lambda & \lambda & \lambda+2G & 0 & 0 & 0\\ 0 & 0 & 0 & \nu & 0 & 0\\ 0 & 0 & 0 & 0 & \nu & 0\\ 0 & 0 & 0 & 0 & 0 & \nu\\ \end{pmatrix}\end{split}\end{aligned}\end{align} \]

where C is stiffness tensor, \(\nu\) is shear modulus, \(\lambda\) is lambda coefficient.

Parameters:
  • deltaEffectiveStress (npt.NDArray[np.float64]) – effective stress variation (\(\Delta\sigma_{eff}\) - Pa) [S11, S22, S33, S23, S13, S12]

  • youngModulus (npt.NDArray[np.float64]) – Young modulus (E - Pa)

  • poissonRatio (npt.NDArray[np.float64]) – Poisson’s ratio (\(\nu\))

Returns:

elastic strain (\(\epsilon\))

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.lambdaCoefficient(youngModulus, poissonRatio)[source]

Compute lambda coefficient from young modulus and Poisson ratio.

\[\lambda = \frac{E*\nu}{(1+\nu)(1-2\nu)}\]
Parameters:
  • youngModulus (npt.NDArray[np.float64]) – Young modulus (E - Pa)

  • poissonRatio (npt.NDArray[np.float64]) – Poisson’s ratio (\(\nu\))

Returns:

lambda coefficient (\(\lambda\))

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.lithostaticStress(depth, density, gravity)[source]

Compute the lithostatic stress.

Parameters:
  • depth (npt.NDArray[np.float64]) – depth from surface - m)

  • density (npt.NDArray[np.float64]) – density of the overburden (kg/m³)

  • gravity (float) – gravity (m²/s)

Returns:

lithostatic stress (Pa)

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.oedometricModulus(Edef, poissonRatio)[source]

Compute Oedometric modulus.

\[M_{oed} = \frac{E_{def}}{1-2\frac{\nu^2}{1-\nu}}\]
Parameters:
  • Edef (npt.NDArray[np.float64]) – Deformation modulus (\(E_{def}\) - Pa)

  • poissonRatio (npt.NDArray[np.float64]) – Poisson’s ratio (\(\nu\))

Returns:

Oedometric modulus (\(M_{oed}\) - Pa)

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.poissonRatio(bulkModulus, shearModulus)[source]

Compute Poisson’s ratio.

\[\nu = \frac{3K-2G}{2(3K+G)}\]
Parameters:
  • bulkModulus (npt.NDArray[np.float64]) – Bulk modulus (K - Pa)

  • shearModulus (npt.NDArray[np.float64]) – Shear modulus (G - Pa)

Returns:

Poisson’s ratio (\(\nu\))

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.reservoirStressPathOed(biotCoefficient, poissonRatio)[source]

Compute reservoir stress path in oedometric conditions.

\[RSP_{oed}=b\frac{1-2\nu}{1-\nu}\]
Parameters:
  • biotCoefficient (npt.NDArray[np.float64]) – biot coefficient (b)

  • poissonRatio (npt.NDArray[np.float64]) – Poisson’s ratio (\(\nu\))

Returns:

reservoir stress path (\(RSP_{oed}\))

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.reservoirStressPathReal(deltaStress, deltaPressure)[source]

Compute real reservoir stress path.

\[RSP_{real}=\frac{\Delta\sigma}{\Delta P}\]
Parameters:
  • deltaStress (npt.NDArray[np.float64]) – stress difference from start (\(\Delta\sigma\) - Pa)

  • deltaPressure (npt.NDArray[np.float64]) – pressure difference from start (\(\Delta P\) - Pa)

Returns:

reservoir stress path (\(RSP_{real}\))

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.shearCapacityUtilization(traction, rockCohesion, frictionAngle)[source]

Compute shear capacity utilization (SCU).

\[SCU = \frac{abs(\tau_1)}{\tau_{max}}\]

where tau_{max} is the Mohr-Coulomb failure threshold.

Parameters:
  • traction (npt.NDArray[np.float64]) – traction vector (\((\sigma, \tau_1, \tau2)\) - Pa)

  • rockCohesion (float) – rock cohesion (c - Pa).

  • frictionAngle (float) – friction angle (\(\alpha\) - rad).

Returns:

SCU

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.shearModulus(youngModulus, poissonRatio)[source]

Compute shear Modulus from young modulus and poisson ratio.

\[G = \frac{E}{2(1+\nu)}\]
Parameters:
  • youngModulus (npt.NDArray[np.float64]) – Young modulus (E - Pa)

  • poissonRatio (npt.NDArray[np.float64]) – Poisson’s ratio (\(\nu\))

Returns:

Shear modulus (G - Pa)

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.specificGravity(density, specificDensity)[source]

Compute the specific gravity.

\[SG = \frac{\rho}{\rho_f}\]
Parameters:
  • density (npt.NDArray[np.float64]) – density (\(\rho\) - kg/m³)

  • specificDensity (float) – fluid density (\(\rho_f\) - kg/m³)

Returns:

specific gravity (SG - no unit)

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.stressRatio(horizontalStress, verticalStress)[source]

Compute horizontal to vertical stress ratio.

\[r = \frac{\sigma_h}{\sigma_v}\]
Parameters:
  • horizontalStress (npt.NDArray[np.float64]) – horizontal stress (\(\sigma_h\) - Pa)

  • verticalStress (npt.NDArray[np.float64]) – vertical stress (\(\sigma_v\) - Pa)

Returns:

stress ratio (\(\sigma\) - Pa)

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.totalStress(effectiveStress, biot, pressure)[source]

Compute total stress from effective stress, pressure, and Biot coeff.

\[\sigma_{tot} = \sigma_{eff}-bP\]
Parameters:
  • effectiveStress (npt.NDArray[np.float64]) – effective stress (\(\sigma_{eff}\) - Pa) using Geos convention

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

  • pressure (npt.NDArray[np.float64]) – Pore pressure (P - Pa)

Returns:

total stress (\(\sigma_{tot}\) - Pa)

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.totalStressRatioThreshold(pressure, horizontalStress)[source]

Compute total stress ratio threshold.

Corresponds to the fracture threshold from Lemgruber-Traby et al (2024). Fracturing can occur in areas where FT > 1. Equals FractureIndex / totalStressRatio. (see Lemgruber-Traby, A., Cacas, M. C., Bonte, D., Rudkiewicz, J. L., Gout, C., & Cornu, T. (2024). Basin modelling workflow applied to the screening of deep aquifers for potential CO2 storage. Geoenergy, geoenergy2024-010. https://doi.org/10.1144/geoenergy2024-010)

\[\sigma_{Th}=\frac{P}{\sigma_h}\]
Parameters:
  • pressure (npt.NDArray[np.float64]) – Pressure (P - Pa)

  • horizontalStress (npt.NDArray[np.float64]) – minimal horizontal total stress (\(\sigma_h\) - Pa) using Geos convention

Returns:

fracture threshold (\(\sigma_{Th}\))

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geomechanicsCalculatorFunctions.youngModulus(bulkModulus, shearModulus)[source]

Compute Young modulus.

\[E = \frac{9K.G}{3K+G}\]
Parameters:
  • bulkModulus (npt.NDArray[np.float64]) – Bulk modulus (K - Pa)

  • shearModulus (npt.NDArray[np.float64]) – Shear modulus (G - Pa)

Returns:

Young modulus (E - Pa)

Return type:

npt.NDArray[np.float64]

geos_posp.processing.geosLogReaderFunctions module

Functions to read and process Geos log.

geos_posp.processing.geosLogReaderFunctions.buildPropertiesNameForComponents(phasesName)[source]

Builds the list of component property names from the list of phases name.

Parameters:

phasesName (list) – [“CO2”,”Water”]

Returns:

[‘Dissolved mass CO2 in CO2’,’Dissolved mass Water in CO2’, ‘Dissolved mass CO2 in Water’,’Dissolved mass Water in Water’]

Return type:

list

geos_posp.processing.geosLogReaderFunctions.buildPropertiesNameForPhases(nameBlock, phasesName)[source]

Replace phase by phase names.

Parameters:
  • nameBlock (str) – “ Mobile phase mass”

  • phasesName (list[str]) – [“CO2”,”Water”]

Returns:

[‘Mobile CO2 mass’, ‘Mobile Water mass’]

Return type:

list[str]

geos_posp.processing.geosLogReaderFunctions.buildPropertiesNameFromGeosProperties(geosProperties, phasesName)[source]

Extracts the property name and its extensions like min, max, average.

Parameters:
  • geosProperties (str) – “ Delta pressure (min, max)”

  • phasesName (list[str]) – [“CO2”,”Water”]

Returns:

[” Delta pressure min”, “ Delta pressure max”]

Return type:

list[str]

geos_posp.processing.geosLogReaderFunctions.buildPropertiesNameNoPhases(nameBlock, extensions='')[source]

From a name block and extensions, builds a list of properties name.

Parameters:
  • nameBlock (str) – “ Delta pressure “

  • extensions (str) – “min, max)”

Returns:

[” Delta pressure min”, “ Delta pressure max”]

Return type:

list

geos_posp.processing.geosLogReaderFunctions.convertValues(propertyNames, propertyValues, propertiesUnit)[source]

Convert properties to the desired units.

Knowing two lists : 1) float numbers that are supposed to be in SI units ; 2) properties name linked to the float numbers. And that these lists are of same dimension, creates a new list of float values converted to a specific unit linked to the property name.

Parameters:
  • propertyNames (list[str]) – list of property names

  • propertyValues (list[float]) – list of property values.

  • propertiesUnit (dict[str, Unit]) – dictionary of desired units for each property {“pressure”: UnitPressure, …, “propertyTypeN”: UnitPropertyTypeN}

Returns:

list of converted values.

Return type:

list[float]

geos_posp.processing.geosLogReaderFunctions.correctZeroValuesInListOfValues(values)[source]

Replace orhphelin 0 values of input list.

If 0 values are found in a list of values, either replace them with the value found before in the list or keep it 0. We suppose that 2 zeros in a row correspond to a continuity of the values, hence keeping it 0.

Parameters:

values (list[float]) – list of ints or floats

Returns:

list of ints or floats

Return type:

list[float]

geos_posp.processing.geosLogReaderFunctions.countNumberLines(filepath)[source]

Reads a file to find the number of lines within it.

Parameters:

filepath (str) – Path to the file.

Returns:

Number of lines in file.

Return type:

int

geos_posp.processing.geosLogReaderFunctions.elementsAreInLog(filepath, elements)[source]

Indicates if input file contains element from input list of string.

To do so, this reads a file and checks at every line if an element was found within the line. If an element is found, it is not checked again. The function returns True only when there is no more element to check.

Parameters:
  • filepath (str) – Path to the file.

  • elements (list[str]) – Every string that needs to be find inside the file.

Return type:

bool

geos_posp.processing.geosLogReaderFunctions.extractAquifer(geosLogLine)[source]

Extracts the name of the aquifer from a Geos log line.

Parameters:

geosLogLine (str) – #expected line : “ TableFunction: aquifer1_pressureInfluence_table”

Raises:

ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

“aquifer1”

Return type:

str

geos_posp.processing.geosLogReaderFunctions.extractFirstIntFromString(string)[source]

Extracts the first int value from a string.

Parameters:

string (str) – A random string.

Returns:

int or None if no int was found.

geos_posp.processing.geosLogReaderFunctions.extractFloatsFromString(line)[source]

Extracts a list of float numbers from a string.

Parameters:

line (str) – A random string.

Returns:

[float1, …, floatN]

Return type:

list[float]

geos_posp.processing.geosLogReaderFunctions.extractLinearIter(geosLogLine)[source]

From a Geos log line, extracts the int value of linear iterations.

Parameters:

geosLogLine (str) – #expected line : “ Linear Solver | Success | Iterations: 23 | Final Rel Res: 5.96636e-05 | Make Restrictor Time: 0 | Compute Auu Time: 0 | SC Filter Time: 0 | Setup Time: 1.5156 s | Solve Time: 0.041093 s”

Raises:
  • ValueError – “Not enough elements to unpack in time line <<” + geosLogLine + “>>.”

  • ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

23

Return type:

int

geos_posp.processing.geosLogReaderFunctions.extractListIntsFromString(string)[source]

Builds a list of int numbers from a string.

Parameters:

string (str) – A random string.

Returns:

[int1, …, intN]

Return type:

list[int]

geos_posp.processing.geosLogReaderFunctions.extractNewtonIter(geosLogLine)[source]

From a Geos log line, extracts the int value of NewtonIter.

Parameters:

geosLogLine (str) – #expected line : “ Attempt: 0, ConfigurationIter: 0, NewtonIter: 0”

Raises:
  • ValueError – “Not enough elements to unpack in time line <<” + geosLogLine + “>>.”

  • ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

NewtonIter

Return type:

int

geos_posp.processing.geosLogReaderFunctions.extractPhaseId(geosLogLine)[source]

Extracts the phase number id from a Geos log line.

Parameters:

geosLogLine (str) – #expected line : “wellControls1: Phase 0 surface volumetric rate: 30.023748128796043 sm3/s”

Raises:
  • ValueError – “Not enough elements to unpack in line <<” + geosLogLine + “>>.”

  • ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

0

Return type:

int

geos_posp.processing.geosLogReaderFunctions.extractPhaseModel(geosLogLine)[source]

Extracts the name of a phase model from a Geos log line.

Parameters:

geosLogLine (str) – #expected line : “ TableFunction: fluid_phaseModel1_PhillipsBrineDensity_table”

Raises:
  • ValueError – “Not enough elements to unpack in line <<” + geosLogLine + “>>.”

  • ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

“PhillipsBrineDensity”

Return type:

str

geos_posp.processing.geosLogReaderFunctions.extractPropertiesFlow(geosLogLine, phasesName)[source]

Extracts flow property from a Geos log line.

Parameters:
  • geosLogLine (str) – #expected line : “compflowStatistics, Reservoir: Delta pressure (min, max): 0, 0 Pa”

  • phasesName (list[str]) – [“CO2”,”Water”]

Raises:
  • ValueError – “Not enough elements to unpack in line <<” + geosLogLine + “>>.”

  • ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

[“Reservoir__DeltaPressureMin”, “Reservoir__DeltaPressureMax”]

Return type:

list[str]

geos_posp.processing.geosLogReaderFunctions.extractPropertiesWell(geosLogLine, wellName, phaseName)[source]

Extracts the well property presented from a Geos log line.

Parameters:
  • geosLogLine (str) – “wellControls1: Phase 0 surface volumetric rate: 30.023748128796043 sm3/s”

  • wellName (str) – well1

  • phaseName (str) – phase1

Returns:

[“Well1_SurfaceVolumetricRatePhase1”]

Return type:

list[str]

geos_posp.processing.geosLogReaderFunctions.extractRegion(geosLogLine)[source]

Extracts the name of the region from a Geos log line.

Parameters:
  • geosLogLine (str) – #expected line : “Adding Object CellElementRegion

  • ObjectManager::Catalog." (named Reservoir from)

Raises:
  • ValueError – “Not enough elements to unpack in line.”

  • ValueError – “An error has occured while parsing a line.”

Returns:

“Reservoir”

Return type:

str

geos_posp.processing.geosLogReaderFunctions.extractStatsName(geosLogLine)[source]

Extracts the name of the computed statistics name from a Geos log line.

Parameters:

geosLogLine (str) – #expected line :”compflowStatistics, Reservoir: Pressure (min, average, max): 2.86419e+07, 2.93341e+07, 3.006e+07 Pa”

Returns:

“compflowStatistics”

Return type:

str

geos_posp.processing.geosLogReaderFunctions.extractTimeAndDt(geosLogLine)[source]

From a Geos log line, extracts the float values of Time and dt.

Parameters:

geosLogLine (str) – #expected lines : “Time: {} years, {} days, {} hrs, {} min, {} s, dt: {} s, Cycle: {}” “Time: {:.2f} years, dt: {} s, Cycle: {}” “Time: {:.2f} days, dt: {} s, Cycle: {}” “Time: {:.2f} hrs, dt: {} s, Cycle: {}” “Time: {:.2f} min, dt: {} s, Cycle: {}” “Time: {:4.2e} s, dt: {} s, Cycle: {}” “Time: {]s, dt:{}s, Cycle: {}”

Raises:

KeyError – “Cannot add time values for tag=<<timeFactor>>”

Returns:

(time, dt)

Return type:

tuple[float]

geos_posp.processing.geosLogReaderFunctions.extractValueAndNameAquifer(geosLogLine)[source]

Extract value and name of the aquifer contained in a Geos log line.

Parameters:

geosLogLine (str) – “FlowSolverBase compositionalMultiphaseFlow (SimuDeck_aquifer_pression_meme.xml, l.28): at time 100s, the <Aquifer> boundary condition ‘aquifer1’ produces a flux of -0.6181975187076816 kg (or moles if useMass=0).”

Returns:

a tuple with the name and the float value. e.g. (“aquifer1”, -0.6181975187076816)

Return type:

tuple[str, float]

geos_posp.processing.geosLogReaderFunctions.extractValuesFlow(geosLogLine)[source]

Extract values from a Geos log line.

Parameters:

geosLogLine (str) – Geos log line “compflowStatistics, Reservoir: ” “Dissolved component mass: { { 0, 0 }, { 0, -6.38235e+10 } } kg”

Returns:

list of values in the line. [0.0, 0.0, 0.0, -6.38235e+10]

Return type:

list[float]

geos_posp.processing.geosLogReaderFunctions.extractValuesWell(geosLogLine, numberProperties)[source]

Extract values from Geos log line and returns them as a list of floats.

The idea here is first to extract all floats values from the line. Now all of them are useful so we need to keep some of them. Luckily, only the last one or two floats are useful. And to determine that, we use the number of well properties found in the line which ranges from one to two.

Parameters:
  • geosLogLine (str) – “Rank 129: well.CO2010: The density of phase 0 at surface conditions is 1.86 kg/sm3.”

  • numberProperties (int) – number of well properties found in the line.

Returns:

value of the property. e.g. [1.86]

Return type:

list[float]

geos_posp.processing.geosLogReaderFunctions.extractWell(geosLogLine)[source]

Extracts the name of the well from a Geos log line.

Parameters:

geosLogLine (str) – #expected line : “ TableFunction: wellControls_ConstantBHP_table” Or “ TableFunction: wellControls_ConstantPhaseRate_table”

Raises:

ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

“wellControls”

Return type:

str

geos_posp.processing.geosLogReaderFunctions.extractWellTags(geosLogLine)[source]

Extracts the list of well property tags available from a Geos log line.

Parameters:

geosLogLine (str) – line from geos log file.

Returns:

list of tags.

Return type:

list[str]

geos_posp.processing.geosLogReaderFunctions.findNumberPhasesSimulation(filepath)[source]

Find the number of phases from Geos log file.

Geos logs do not have explicit message telling you how many phases were used to perform the simulation, unlike regions, wells etc … Therefore, we need at least to identify the exact number of phases that can be find in this Geos log file to extract correctly properties regarding phase related data.

Parameters:

filepath (str) – Filepath to a Geos log file.

Returns:

The number of phases found in the Geos log.

Return type:

int

geos_posp.processing.geosLogReaderFunctions.formatPropertyName(propertyName)[source]

Clean the string by replacing special characters and removing spaces.

Parameters:

propertyName (str) – name;of:the property

Returns:

NameOfTheProperty

Return type:

str

geos_posp.processing.geosLogReaderFunctions.identifyCurrentWell(geosLogLine, lastWellName)[source]

Identify the current name of the well rom a Geos log line.

Because properties values of wells can be output without specifying the name of the well which they belong to, we have to assume that the name of well for the current properties observed is either :

  • the name of the well that can be found inside the line

  • the name of the last well name found in former lines.

Parameters:
  • geosLogLine (str) – line from Geos log file #expected lines with well name : “Rank 18: well.CO2001: BHP (at the specified reference elevation): 19318538.400682557 Pa” Or “wellControls1: BHP (at the specified reference elevation): 12337146.157562563 Pa” #line with no well name : “The total rate is 0 kg/s, which corresponds to a total surface volumetric rate of 0 sm3/s”

  • lastWellName (str) – name of the last well found

Raises:

ValueError – “An error has occured while parsing <<geosLogLine>>.”

Returns:

“wellControls”

Return type:

str

geos_posp.processing.geosLogReaderFunctions.identifyProperties(properties)[source]

Identify properties and add identifer.

From a list of properties name, identifies each of them with a certain integer, to link it to a meta property by adding an id in front of the property name.

Parameters:

properties (list[str]) – [“CaprockPressureMax”, “CaprockPressureMin”]

Returns:

[1:”CaprockPressureMax”, 1:”CaprockPressureMin”]

Return type:

list[tuple[str, int]]

geos_posp.processing.geosLogReaderFunctions.isFloat(element)[source]

Check whether an element is float or not.

Parameters:

element (Any) – input number to test.

Returns:

True if the number is a float.

Return type:

bool

geos_posp.processing.geosLogReaderFunctions.phaseNamesBuilder(numberPhases, phasesFromUser)[source]

Build phase names.

When creating phase names, the user can or cannot have defined his own phase names when reading the log. Therefore, whether phase names were provided or not, if we have N phases found in the log, a list of N phase names will be created, starting from phase0 up to phaseN.

Parameters:
  • numberPhases (int) – Number of phases found in the log file.

  • phasesFromUser (list[str]) – names chosen by the user, can be more or less than numberPhases.

Returns:

[nameFromUser0, nameFromUser1, …, phaseN-1, phaseN]

Return type:

list[str]

geos_posp.processing.geosLogReaderFunctions.replaceSpecialCharactersWithWhitespace(sentence)[source]

Replace every special characters in a string with whitespaces.

Parameters:

sentence (str) – Random string “hi ‘(_there(‘’&*$^,:;’”

Returns:

“hi there “

Return type:

str

geos_posp.processing.geosLogReaderFunctions.timeInSecond(timeCounter)[source]

Calculates the time in s from a dict of different time quantities.

Parameters:

timeCounter (dict[str, float]) – timeCounter: {“years”: x0, “days”: x0, “hrs”: x0, “min”: x0, “s”: x0}

Returns:

Sum in seconds of all time quantities.

Return type:

float

geos_posp.processing.geosLogReaderFunctions.transformUserChoiceToListPhases(userChoice)[source]

Get a list of phase name from the input string.

When using GeosLogReader, the user can choose the names of the phases to use. The wished format is to specify each name in the good, separated by whitespaces or either commas.

Parameters:

userChoice (str | None) – Output from EnterPhaseNames string vector widget.

Returns:

[phase0, phase1, …, phaseN]

Return type:

list[str]

geos_posp.processing.MohrCircle module

MohrCircle module define the Mohr’s circle parameters.

Inputs are a 6 component stress vector, a circle id, and the mechanical convention used for compression. The class computes principal components from stress vector during initialization. Accessors get access to these 3 principal components as well as circle center and radius.

To use the object:

from processing.MohrCircle import MohrCircle

# create the object
stressVector :npt.NDArray[np.float64]
circleId :str
mohrCircle :MohrCircle = MohrCircle(circleId)

# either directly set principal components (p3 <= p2 <= p1)
mohrCircle.SetPrincipalComponents(p3, p2, p1)
# or compute them from stress vector
mohrCircle.computePrincipalComponents(stressVector)

# access to members
id :str = mohrCircle.getCircleId()
p1, p2, p3 :float = mohrCircle.getPrincipalComponents()
radius :float = mohrCircle.getCircleRadius()
center :float = mohrCircle.getCircleCenter()
class geos_posp.processing.MohrCircle.MohrCircle(circleId)[source]

Bases: object

Compute Mohr’s Circle from input stress.

Parameters:

circleId (str) – Mohr’s circle id.

computePrincipalComponents(stressVector)[source]

Calculate principal components.

Parameters:

stressVector (npt.NDArray[np.float64]) – 6 components stress vector Stress vector must follow GEOS convention (XX, YY, ZZ, YZ, XZ, XY)

getCircleCenter()[source]

Compute and return Mohr’s circle center from principal components.

getCircleId()[source]

Access the Id of the Mohr circle.

Returns:

Id of the Mohr circle

Return type:

str

getCircleRadius()[source]

Compute and return Mohr’s circle radius from principal components.

getPrincipalComponents()[source]

Get Moh’s circle principal components.

setCircleId(circleId)[source]

Set circle Id variable.

Parameters:

circleId (str) – circle Id.

setPrincipalComponents(p3, p2, p1)[source]

Set principal components.

Parameters:
  • p3 (float) – first component. Must be the lowest.

  • p2 (float) – second component.

  • p1 (float) – third component. Must be the greatest.

geos_posp.processing.MohrCoulomb module

MohrCoulomb module define the Mohr-Coulomb failure envelop class.

Inputs are the rock cohesion (Pa) and the friction angle (°). 2 methods allow to compute either shear stress values according to normal stress values, or the failure envelope including normal stress and corresponding shear stress values.

To use the object:

from processing.MohrCoulomb import MohrCoulomb

# create the object
rockCohesion :float = 1.5e9 # Pa
frictionAngle :float = 10 # degree
mohrCoulomb = MohrCoulomb(rockCohesion, frictionAngle)

# compute shear stress values according to normal stress values
normalStress :npt.NDArray[np.float64] = np.linspace(1e9, 1.5e9)
shearStress :npt.NDArray[np.float64] = mohrCoulomb.computeShearStress(normalStress)

# compute the failure envelope including normal stress and corresponding shear
# stress values
# ones may also define minimum normal stress and the number of points
normalStressMax :float = 1.5e9
normalStress, shearStress = mohrCoulomb.computeFailureEnvelop(normalStressMax)
class geos_posp.processing.MohrCoulomb.MohrCoulomb(rockCohesion, frictionAngle)[source]

Bases: object

Define Mohr-Coulomb failure envelop.

Parameters:
  • rockCohesion (float) – rock cohesion (Pa).

  • frictionAngle (float) – friction angle (rad).

computeFailureEnvelop(stressNormalMax, stressNormalMin=None, n=100)[source]

Compute the envelope failure between min and max normal stress.

Parameters:
  • stressNormalMax (float) – Maximum normal stress (Pa)

  • stressNormalMin (float | None, optional) –

    Minimum normal stress. If it is None, the envelop is computed from the its intersection with the abscissa axis.

    Defaults to None.

  • n (int, optional) – Number of points to define the envelop.

  • 100. (Defaults to)

Returns:

failure envelop coordinates where first element is the abscissa and second element is the ordinates.

Return type:

tuple[npt.NDArray[np.float64], npt.NDArray[np.float64]]

computeShearStress(stressNormal0)[source]

Compute shear stress from normal stress.

Parameters:

stressNormal0 (float | npt.NDArray[np.float64]) – normal stress value (Pa)

Returns:

shear stress value.

Return type:

float | npt.NDArray[np.float64])

geos_posp.processing.multiblockInpectorTreeFunctions module

Functions to explore and process multiblock inspector tree.

geos_posp.processing.multiblockInpectorTreeFunctions.getBlockElementIndexes(input)[source]

Get a list of list that contains flat indexes of elementary blocks.

Each sublist contains the indexes of elementary blocks that belongs to a same parent node.

Parameters:

input (vtkMultiBlockDataSet | vtkCompositeDataSet) – input multi block object.

Returns:

list of list of flat indexes

Return type:

list[list[int]]

geos_posp.processing.multiblockInpectorTreeFunctions.getBlockElementIndexesFlatten(input)[source]

Get a flatten list that contains flat indexes of elementary blocks.

Parameters:

input (vtkMultiBlockDataSet | vtkCompositeDataSet) – input multi block object.

Returns:

list of flat indexes

Return type:

list[int]

geos_posp.processing.multiblockInpectorTreeFunctions.getBlockFromFlatIndex(multiBlock, blockIndex)[source]

Get the block with blockIndex from input vtkMultiBlockDataSet.

Parameters:
  • multiBlock (vtkMultiBlockDataSet | vtkCompositeDataSet) – input multi block

  • blockIndex (int) – block index

Returns:

block if it exists, None otherwise

Return type:

vtkMultiBlockDataSet

geos_posp.processing.multiblockInpectorTreeFunctions.getBlockFromName(multiBlock, blockName)[source]

Get the block named blockName from input vtkMultiBlockDataSet.

Parameters:
  • multiBlock (vtkMultiBlockDataSet | vtkCompositeDataSet) – input multi block

  • blockName (str) – block name

Returns:

block if it exists, None otherwise

Return type:

vtkDataObject

geos_posp.processing.multiblockInpectorTreeFunctions.getBlockIndexFromName(input, name)[source]

Get block flat index from name of node in the vtkMultiBlockDataSet tree.

Parameters:
  • input (vtkMultiBlockDataSet | vtkCompositeDataSet) – input multi block object.

  • name (str) – name of the block to get the index in the tree.

Returns:

index of the block if found, -1 otherwise.

Return type:

int

geos_posp.processing.multiblockInpectorTreeFunctions.getBlockName(input)[source]

Get the name of input block.

If input is a vtkMultiBlockDataSet or vtkCompositeDataSet, returns the name of the lowest level unique child block.

Parameters:

input (vtkMultiBlockDataSet | vtkCompositeDataSet) – input multi block object.

Returns:

name of the block in the tree.

Return type:

str

geos_posp.processing.multiblockInpectorTreeFunctions.getBlockNameFromIndex(input, index)[source]

Get the name of a block from input multiblock and block flat index.

Parameters:
  • input (vtkMultiBlockDataSet | vtkCompositeDataSet) – input multi block object.

  • index (int) – flat index of the block to get the name

Returns:

name of the block in the tree.

Return type:

str

geos_posp.processing.multiblockInpectorTreeFunctions.getElementaryCompositeBlockIndexes(input)[source]

Get indexes of composite block that contains elementrary blocks.

Parameters:

input (vtkMultiBlockDataSet | vtkCompositeDataSet) – input multi block object.

Returns:

dictionary that contains names as keys and flat indices as values of the parent composite blocks of elementary blocks.

Return type:

dict[str, int]

geos_posp.processing.vtkUtils module

Utilities to process vtk objects.

geos_posp.processing.vtkUtils.computeCellCenterCoordinates(mesh)[source]

Get the coordinates of Cell center.

Parameters:

mesh (vtkDataSet) – input surface

Returns:

cell center coordinates

Return type:

vtkPoints

geos_posp.processing.vtkUtils.copyAttribute(objectFrom, objectTo, attributNameFrom, attributNameTo)[source]

Copy an attribute from objectFrom to objectTo.

Parameters:
  • objectFrom (vtkMultiBlockDataSet) – object from which to copy the attribute.

  • objectTo (vtkMultiBlockDataSet) – object where to copy the attribute.

  • attributNameFrom (str) – attribute name in objectFrom.

  • attributNameTo (str) – attribute name in objectTo.

Returns:

True if copy sussfully ended, False otherwise

Return type:

bool

geos_posp.processing.vtkUtils.copyAttributeDataSet(objectFrom, objectTo, attributNameFrom, attributNameTo)[source]

Copy an attribute from objectFrom to objectTo.

Parameters:
  • objectFrom (vtkDataSet) – object from which to copy the attribute.

  • objectTo (vtkDataSet) – object where to copy the attribute.

  • attributNameFrom (str) – attribute name in objectFrom.

  • attributNameTo (str) – attribute name in objectTo.

Returns:

True if copy sussfully ended, False otherwise

Return type:

bool

geos_posp.processing.vtkUtils.createAttribute(dataSet, array, attributeName, componentNames, onPoints)[source]

Create an attribute from the given array.

Parameters:
  • dataSet (vtkDataSet) – dataSet where to create the attribute

  • array (npt.NDArray[np.float64]) – array that contains the values

  • attributeName (str) – name of the attribute

  • componentNames (tuple[str,...]) – name of the components for vectorial attributes

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

True if the attribute was correctly created

Return type:

bool

geos_posp.processing.vtkUtils.createCellCenterAttribute(mesh, cellCenterAttributeName)[source]

Create elementCenter attribute if it does not exist.

Parameters:
  • mesh (vtkMultiBlockDataSet | vtkDataSet) – input mesh

  • cellCenterAttributeName (str) – Name of the attribute

Raises:
  • TypeError – Raised if input mesh is not a vtkMultiBlockDataSet or a

  • vtkDataSet.

Returns:

True if calculation successfully ended, False otherwise.

Return type:

bool

geos_posp.processing.vtkUtils.createConstantAttribute(object, values, attributeName, componentNames, onPoints)[source]

Create an attribute with a constant value everywhere if absent.

Parameters:
  • object (vtkDataObject) – object (vtkMultiBlockDataSet, vtkDataSet) where to create the attribute

  • values (list[float]) – list of values of the attribute for each components

  • attributeName (str) – name of the attribute

  • componentNames (tuple[str,...]) – name of the components for vectorial attributes

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

True if the attribute was correctly created

Return type:

bool

geos_posp.processing.vtkUtils.createConstantAttributeDataSet(dataSet, values, attributeName, componentNames, onPoints)[source]

Create an attribute with a constant value everywhere.

Parameters:
  • dataSet (vtkDataSet) – vtkDataSet where to create the attribute

  • values (list[float]) – list of values of the attribute for each components

  • attributeName (str) – name of the attribute

  • componentNames (tuple[str,...]) – name of the components for vectorial attributes

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

True if the attribute was correctly created

Return type:

bool

geos_posp.processing.vtkUtils.createConstantAttributeMultiBlock(multiBlockDataSet, values, attributeName, componentNames, onPoints)[source]

Create an attribute with a constant value everywhere if absent.

Parameters:
  • multiBlockDataSet (vtkMultiBlockDataSet | vtkCompositeDataSet) – vtkMultiBlockDataSet where to create the attribute

  • values (list[float]) – list of values of the attribute for each components

  • attributeName (str) – name of the attribute

  • componentNames (tuple[str,...]) – name of the components for vectorial attributes

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

True if the attribute was correctly created

Return type:

bool

geos_posp.processing.vtkUtils.createEmptyAttribute(object, attributeName, componentNames, dataType, onPoints)[source]

Create an empty attribute.

Parameters:
  • object (vtkDataObject) – object (vtkMultiBlockDataSet, vtkDataSet) where to create the attribute

  • attributeName (str) – name of the attribute

  • componentNames (tuple[str,...]) – name of the components for vectorial attributes

  • dataType (int) – data type.

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

True if the attribute was correctly created

Return type:

bool

geos_posp.processing.vtkUtils.doCreateCellCenterAttribute(block, cellCenterAttributeName)[source]

Create elementCenter attribute in a vtkDataSet if it does not exist.

Parameters:
  • block (vtkDataSet) – input mesh that must be a vtkDataSet

  • cellCenterAttributeName (str) – Name of the attribute

Returns:

True if calculation successfully ended, False otherwise.

Return type:

bool

geos_posp.processing.vtkUtils.extractBlock(multiBlockDataSet, blockIndex)[source]

Extract the block with index blockIndex from multiBlockDataSet.

Parameters:
  • multiBlockDataSet (vtkMultiBlockDataSet) – multiblock dataset from which to extract the block

  • blockIndex (int) – block index to extract

Returns:

extracted block

Return type:

vtkMultiBlockDataSet

geos_posp.processing.vtkUtils.extractSurfaceFromElevation(mesh, elevation)[source]

Extract surface at a constant elevation from a mesh.

Parameters:
  • mesh (vtkUnstructuredGrid) – input mesh

  • elevation (float) – elevation at which to extract the surface

Returns:

output surface

Return type:

vtkPolyData

geos_posp.processing.vtkUtils.fillAllPartialAttributes(multiBlockMesh, onPoints=False)[source]

Fill all the partial attributes of multiBlockMesh with nan values.

Parameters:
  • multiBlockMesh (vtkMultiBlockDataSet | vtkCompositeDataSet | vtkDataObject) – multiBlockMesh where to fill the attribute

  • onPoints (bool, optional) –

    Attribute is on Points (False) or on Cells.

    Defaults to False.

Returns:

True if calculation successfully ended, False otherwise

Return type:

bool

geos_posp.processing.vtkUtils.fillPartialAttributes(multiBlockMesh, attributeName, nbComponents, onPoints=False)[source]

Fill input partial attribute of multiBlockMesh with nan values.

Parameters:
  • multiBlockMesh (vtkMultiBlockDataSet | vtkCompositeDataSet | vtkDataObject) – multiBlock mesh where to fill the attribute

  • attributeName (str) – attribute name

  • nbComponents (int) – number of components

  • onPoints (bool, optional) –

    Attribute is on Points (False) or on Cells.

    Defaults to False.

Returns:

True if calculation successfully ended, False otherwise

Return type:

bool

geos_posp.processing.vtkUtils.getArrayInObject(object, attributeName, onPoints)[source]

Return the numpy array corresponding to input attribute name in table.

Parameters:
  • object (PointSet or UnstructuredGrid) – input object

  • attributeName (str) – name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

the array corresponding to input attribute name.

Return type:

ArrayLike[float]

geos_posp.processing.vtkUtils.getAttributeSet(object, onPoints)[source]

Get the set of all attributes from an object on points or on cells.

Parameters:
  • object (Any) – object where to find the attributes.

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

set of attribute names present in input object.

Return type:

set[str]

geos_posp.processing.vtkUtils.getAttributeValuesAsDF(surface, attributeNames)[source]

Get attribute values from input surface.

Parameters:
  • surface (vtkPolyData) – mesh where to get attribute values

  • attributeNames (tuple[str,...]) – tuple of attribute names to get the values.

Returns:

DataFrame containing property names as columns.

Return type:

pd.DataFrame

geos_posp.processing.vtkUtils.getAttributesFromDataSet(object, onPoints)[source]

Get the dictionnary of all attributes of a vtkDataSet on points or cells.

Parameters:
  • object (vtkDataSet) – object where to find the attributes.

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

List of the names of the attributes.

Return type:

dict[str, int]

geos_posp.processing.vtkUtils.getAttributesFromMultiBlockDataSet(object, onPoints)[source]

Get the dictionnary of all attributes of object on points or on cells.

Parameters:
  • object (vtkMultiBlockDataSet | vtkCompositeDataSet) – object where to find the attributes.

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

Dictionnary of the names of the attributes as keys, and

number of components as values.

Return type:

dict[str, int]

geos_posp.processing.vtkUtils.getAttributesWithNumberOfComponents(object, onPoints)[source]

Get the dictionnary of all attributes from object on points or cells.

Parameters:
  • object (Any) – object where to find the attributes.

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

dictionnary where keys are the names of the attributes

and values the number of components.

Return type:

dict[str, int]

geos_posp.processing.vtkUtils.getBounds(input)[source]

Get bounds of either single of composite data set.

Parameters:

input (Union[vtkUnstructuredGrid, vtkMultiBlockDataSet]) – input mesh

Returns:

tuple containing

bounds (xmin, xmax, ymin, ymax, zmin, zmax)

Return type:

tuple[float, float, float, float, float, float]

geos_posp.processing.vtkUtils.getComponentNames(dataSet, attributeName, onPoints)[source]

Get the name of the components of attribute attributeName in dataSet.

Parameters:
  • dataSet (vtkDataSet | vtkMultiBlockDataSet | vtkCompositeDataSet | vtkDataObject) – dataSet where the attribute is.

  • attributeName (str) – name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

names of the components.

Return type:

tuple[str,…]

geos_posp.processing.vtkUtils.getComponentNamesDataSet(dataSet, attributeName, onPoints)[source]

Get the name of the components of attribute attributeName in dataSet.

Parameters:
  • dataSet (vtkDataSet) – dataSet where the attribute is.

  • attributeName (str) – name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

names of the components.

Return type:

tuple[str,…]

geos_posp.processing.vtkUtils.getComponentNamesMultiBlock(dataSet, attributeName, onPoints)[source]

Get the name of the components of attribute in MultiBlockDataSet.

Parameters:
  • dataSet (vtkMultiBlockDataSet | vtkCompositeDataSet) – dataSet where the attribute is.

  • attributeName (str) – name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

names of the components.

Return type:

tuple[str,…]

geos_posp.processing.vtkUtils.getMonoBlockBounds(input)[source]

Get boundary box extrema coordinates for a vtkUnstructuredGrid.

Parameters:

input (vtkMultiBlockDataSet) – input single block mesh

Returns:

tuple containing

bounds (xmin, xmax, ymin, ymax, zmin, zmax)

Return type:

tuple[float, float, float, float, float, float]

geos_posp.processing.vtkUtils.getMultiBlockBounds(input)[source]

Get boundary box extrema coordinates for a vtkMultiBlockDataSet.

Parameters:

input (vtkMultiBlockDataSet) – input multiblock mesh

Returns:

bounds.

Return type:

tuple[float, float, float, float, float, float]

geos_posp.processing.vtkUtils.getNumberOfComponents(dataSet, attributeName, onPoints)[source]

Get the number of components of attribute attributeName in dataSet.

Parameters:
  • dataSet (vtkMultiBlockDataSet | vtkCompositeDataSet | vtkDataSet) – dataSet where the attribute is.

  • attributeName (str) – name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

number of components.

Return type:

int

geos_posp.processing.vtkUtils.getNumberOfComponentsDataSet(dataSet, attributeName, onPoints)[source]

Get the number of components of attribute attributeName in dataSet.

Parameters:
  • dataSet (vtkDataSet) – dataSet where the attribute is.

  • attributeName (str) – name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

number of components.

Return type:

int

geos_posp.processing.vtkUtils.getNumberOfComponentsMultiBlock(dataSet, attributeName, onPoints)[source]

Get the number of components of attribute attributeName in dataSet.

Parameters:
  • dataSet (vtkMultiBlockDataSet | vtkCompositeDataSet) – multi block data Set where the attribute is.

  • attributeName (str) – name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

number of components.

Return type:

int

geos_posp.processing.vtkUtils.getVtkArrayInObject(object, attributeName, onPoints)[source]

Return the array corresponding to input attribute name in table.

Parameters:
  • object (PointSet or UnstructuredGrid) – input object

  • attributeName (str) – name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

the vtk array corresponding to input attribute name.

Return type:

vtkDoubleArray

geos_posp.processing.vtkUtils.isAttributeInObject(object, attributeName, onPoints)[source]

Check if an attribute is in the input object.

Parameters:
  • object (vtkMultiBlockDataSet | vtkDataSet) – input object

  • attributeName (str) – name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

True if the attribute is in the table, False otherwise

Return type:

bool

geos_posp.processing.vtkUtils.isAttributeInObjectDataSet(object, attributeName, onPoints)[source]

Check if an attribute is in the input object.

Parameters:
  • object (vtkDataSet) – input object

  • attributeName (str) – name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

True if the attribute is in the table, False otherwise

Return type:

bool

geos_posp.processing.vtkUtils.isAttributeInObjectMultiBlockDataSet(object, attributeName, onPoints)[source]

Check if an attribute is in the input object.

Parameters:
  • object (vtkMultiBlockDataSet) – input multiblock object

  • attributeName (str) – name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

True if the attribute is in the table, False otherwise

Return type:

bool

geos_posp.processing.vtkUtils.mergeBlocks(input, keepPartialAttributes=False)[source]

Merge all blocks of a multi block mesh.

Parameters:
  • input (vtkMultiBlockDataSet | vtkCompositeDataSet) – composite object to merge blocks

  • keepPartialAttributes (bool) –

    if True, keep partial attributes after merge.

    Defaults to False.

Returns:

merged block object

Return type:

vtkUnstructuredGrid

geos_posp.processing.vtkUtils.renameAttribute(object, attributeName, newAttributeName, onPoints)[source]

Rename an attribute.

Parameters:
  • object (vtkMultiBlockDataSet) – object where the attribute is

  • attributeName (str) – name of the attribute

  • newAttributeName (str) – new name of the attribute

  • onPoints (bool) – True if attributes are on points, False if they are on cells.

Returns:

True if renaming operation successfully ended.

Return type:

bool

geos_posp.processing.vtkUtils.transferPointDataToCellData(mesh)[source]

Transfer point data to cell data.

Parameters:

mesh (vtkPointSet) – Input mesh.

Returns:

Output mesh where point data were transferred to cells.

Return type:

vtkPointSet