PyHydroGeophysX.petrophysics package

Submodules

PyHydroGeophysX.petrophysics.resistivity_models module

Simplified Waxman-Smits model for converting between water content and resistivity.

This implementation follows the Waxman-Smits model that expresses conductivity as:

σ = σsat * S^n + σs * S^(n-1)

where: - σ is the electrical conductivity of the formation - σsat is the conductivity at full saturation without surface effects (1/rhos) - σs is the surface conductivity - S is the water saturation (S = θ/φ where θ is water content and φ is porosity) - n is the saturation exponent

The resistivity is the reciprocal of conductivity: ρ = 1/σ

PyHydroGeophysX.petrophysics.resistivity_models.WS_Model(saturation, porosity, sigma_w, m, n, sigma_s=0)[source]

Convert water content to resistivity using Waxman-Smits model.

Based on equation: σ = (S_w^n/F)σ_w + σ_s where F = φ^(-m) (formation factor from Archie’s law)

Parameters:
  • saturation (array) – Saturation (S_w)

  • porosity (array) – Porosity values (φ)

  • sigma_w (float) – Pore water conductivity (σ_w)

  • m (float) – Cementation exponent

  • n (float) – Saturation exponent

  • sigma_s (float) – Surface conductivity (σ_s). Default is 0 (no surface effects).

Returns:

Resistivity values

Return type:

array

PyHydroGeophysX.petrophysics.resistivity_models.resistivity_to_porosity(resistivity, saturation, m, rho_fluid, n, sigma_sur=0, a=1.0)[source]

Convert resistivity to porosity using Waxman-Smits model, given known saturation.

The function solves the Waxman-Smits equation for porosity: 1/rho = sigma_sat * S^n + sigma_sur * S^(n-1) where sigma_sat = 1/rhos and rhos = a * rho_fluid * porosity^(-m)

Rearranging: porosity = [(1/rho - sigma_sur * S^(n-1)) * (a * rho_fluid) / S^n]^(1/m)

Parameters:
  • resistivity (array) – Resistivity values (ohm-m)

  • saturation (array) – Saturation values (fraction, 0-1)

  • m (float) – Cementation exponent (typically 1.3-2.5)

  • rho_fluid (float) – Fluid resistivity (ohm-m)

  • n (float) – Saturation exponent (typically 1.8-2.2)

  • sigma_sur (float) – Surface conductivity (S/m). Default is 0 (no surface effects)

  • a (float) – Tortuosity factor. Default is 1.0

Returns:

Porosity values (fraction, 0-1)

Return type:

array

PyHydroGeophysX.petrophysics.resistivity_models.resistivity_to_saturation(resistivity, porosity, m, rho_fluid, n, sigma_sur=0, a=1.0)[source]

Convert resistivity to saturation using Waxman-Smits model.

The function calculates saturated resistivity using Archie’s law: rhos = a * rho_fluid * porosity^(-m)

Then solves the Waxman-Smits equation: 1/rho = sigma_sat * S^n + sigma_sur * S^(n-1) where sigma_sat = 1/rhos

Parameters:
  • resistivity (array) – Resistivity values (ohm-m)

  • porosity (array) – Porosity values (fraction, 0-1)

  • m (float) – Cementation exponent (typically 1.3-2.5)

  • rho_fluid (float) – Fluid resistivity (ohm-m)

  • n (float) – Saturation exponent (typically 1.8-2.2)

  • sigma_sur (float) – Surface conductivity (S/m). Default is 0 (no surface effects)

  • a (float) – Tortuosity factor. Default is 1.0

Returns:

Saturation values (fraction, 0-1)

Return type:

array

PyHydroGeophysX.petrophysics.resistivity_models.resistivity_to_saturation2(resistivity, rhos, n, sigma_sur=0)[source]

Convert resistivity to saturation using Waxman-Smits model.

Parameters:
  • resistivity (array) – Resistivity values

  • rhos (float) – Saturated resistivity without surface effects

  • n (float) – Saturation exponent

  • sigma_sur (float) – Surface conductivity. Default is 0 (no surface effects).

Returns:

Saturation values

Return type:

array

PyHydroGeophysX.petrophysics.resistivity_models.resistivity_to_water_content(resistivity, rhos, n, porosity, sigma_sur=0)[source]

Convert resistivity to water content using Waxman-Smits model.

Parameters:
  • resistivity (array) – Resistivity values

  • rhos (float) – Saturated resistivity without surface effects

  • n (float) – Saturation exponent

  • porosity (array) – Porosity values

  • sigma_sur (float) – Surface conductivity. Default is 0 (no surface effects).

Returns:

Volumetric water content values

Return type:

array

PyHydroGeophysX.petrophysics.resistivity_models.water_content_to_resistivity(water_content, rhos, n, porosity, sigma_sur=0)[source]

Convert water content to resistivity using Waxman-Smits model.

Parameters:
  • water_content (array) – Volumetric water content (θ)

  • rhos (float) – Saturated resistivity without surface effects

  • n (float) – Saturation exponent

  • porosity (array) – Porosity values (φ)

  • sigma_sur (float) – Surface conductivity. Default is 0 (no surface effects).

Returns:

Resistivity values

Return type:

array

PyHydroGeophysX.petrophysics.velocity_models module

Seismic velocity models for relating rock properties to elastic wave velocities.

class PyHydroGeophysX.petrophysics.velocity_models.BaseVelocityModel[source]

Bases: object

Base class for seismic velocity models.

calculate_velocity(**kwargs) ndarray[source]

Calculate seismic velocity from rock properties.

Parameters:

**kwargs – Rock properties specific to each model

Returns:

Seismic velocity values (Vp, Vs, or both)

class PyHydroGeophysX.petrophysics.velocity_models.BrieModel(exponent: float = 3.0)[source]

Bases: object

Brie’s model for calculating the effective bulk modulus of a partially saturated medium.

calculate_fluid_modulus(saturation: float, water_modulus: float = 2.0, gas_modulus: float = 0.01) float[source]

Calculate effective fluid bulk modulus using Brie’s equation.

Parameters:
  • saturation – Water saturation (0 to 1)

  • water_modulus – Bulk modulus of water (GPa, default: 2.0)

  • gas_modulus – Bulk modulus of gas (GPa, default: 0.01)

Returns:

Effective fluid bulk modulus (GPa)

calculate_saturated_modulus(dry_modulus: float, mineral_modulus: float, porosity: float, saturation: float, water_modulus: float = 2.0, gas_modulus: float = 0.01) float[source]

Calculate the saturated bulk modulus based on Brie’s equation.

Parameters:
  • dry_modulus – Bulk modulus of the dry rock (GPa)

  • mineral_modulus – Bulk modulus of the mineral matrix (GPa)

  • porosity – Porosity of the rock

  • saturation – Water saturation (0 to 1)

  • water_modulus – Bulk modulus of water (GPa, default: 2.0)

  • gas_modulus – Bulk modulus of gas (GPa, default: 0.01)

Returns:

Saturated bulk modulus (GPa)

class PyHydroGeophysX.petrophysics.velocity_models.DEMModel[source]

Bases: BaseVelocityModel

Differential Effective Medium (DEM) model for calculating elastic properties and seismic velocities of porous rocks.

calculate_velocity(porosity: ndarray, saturation: ndarray, bulk_modulus: float, shear_modulus: float, mineral_density: float, aspect_ratio: float = 0.1) Tuple[ndarray, ndarray, ndarray][source]

Calculate P-wave velocity using the DEM model.

Parameters:
  • porosity – Porosity values (array)

  • saturation – Saturation values (array)

  • bulk_modulus – Initial bulk modulus of the solid matrix (GPa)

  • shear_modulus – Initial shear modulus of the solid matrix (GPa)

  • mineral_density – Density of the solid matrix (kg/m³)

  • aspect_ratio – Aspect ratio of pores (default: 0.1)

Returns:

Effective bulk modulus (GPa), effective shear modulus (GPa), and P-wave velocity (m/s)

class PyHydroGeophysX.petrophysics.velocity_models.HertzMindlinModel(critical_porosity: float = 0.4, coordination_number: float = 4.0)[source]

Bases: BaseVelocityModel

Hertz-Mindlin model and Hashin-Shtrikman bounds for seismic velocity in porous rocks.

calculate_velocity(porosity: ndarray, saturation: ndarray, bulk_modulus: float, shear_modulus: float, mineral_density: float, depth: float = 1.0) Tuple[ndarray, ndarray][source]

Calculate P-wave velocity for porous rocks.

Parameters:
  • porosity – Porosity values (array)

  • saturation – Saturation values (array)

  • bulk_modulus – Bulk modulus of the solid matrix (GPa)

  • shear_modulus – Shear modulus of the solid matrix (GPa)

  • mineral_density – Density of the solid matrix (kg/m³)

  • depth – Depth for pressure estimation (m, default: 1.0)

Returns:

Tuple of high-bound P-wave velocity (m/s) and low-bound P-wave velocity (m/s)

class PyHydroGeophysX.petrophysics.velocity_models.VRHModel[source]

Bases: BaseVelocityModel

Voigt-Reuss-Hill (VRH) mixing model for effective elastic properties of composites.

calculate_properties(fractions: List[float], bulk_moduli: List[float], shear_moduli: List[float], densities: List[float]) Tuple[float, float, float][source]

Calculate effective elastic properties using the VRH model.

Parameters:
  • fractions – Volume fractions of each mineral (must sum to 1)

  • bulk_moduli – Bulk moduli of each mineral (GPa)

  • shear_moduli – Shear moduli of each mineral (GPa)

  • densities – Densities of each mineral (kg/m³)

Returns:

Effective bulk modulus (GPa), effective shear modulus (GPa), and effective density (kg/m³)

calculate_velocity(fractions: List[float], bulk_moduli: List[float], shear_moduli: List[float], densities: List[float]) Tuple[float, float][source]

Calculate P-wave and S-wave velocities using the VRH model.

Parameters:
  • fractions – Volume fractions of each mineral (must sum to 1)

  • bulk_moduli – Bulk moduli of each mineral (GPa)

  • shear_moduli – Shear moduli of each mineral (GPa)

  • densities – Densities of each mineral (kg/m³)

Returns:

P-wave velocity (m/s) and S-wave velocity (m/s)

PyHydroGeophysX.petrophysics.velocity_models.VRH_model(f=[0.35, 0.25, 0.2, 0.125, 0.075], K=[55.4, 36.6, 75.6, 46.7, 50.4], G=[28.1, 45, 25.6, 23.65, 27.4], rho=[2560, 2650, 2630, 2540, 3050])[source]

Implements the Voigt-Reuss-Hill (VRH) mixing model to estimate the effective bulk modulus (Km), shear modulus (Gm), and density (rho_b) of a composite material made from various minerals.

Parameters: f (list): Fraction of each mineral in the composite (must sum to 1). K (list): Bulk modulus of each mineral (GPa). G (list): Shear modulus of each mineral (GPa). rho (list): Density of each mineral (kg/m^3).

Returns: Km (float): Effective bulk modulus of the composite material (GPa). Gm (float): Effective shear modulus of the composite material (GPa). rho_b (float): Effective density of the composite material (kg/m^3).

PyHydroGeophysX.petrophysics.velocity_models.satK(Keff, Km, phi, Sat)[source]

Calculate the saturated bulk modulus (K_sat) based on Brie’s equation.

Parameters: Keff (float): Effective bulk modulus of the dry rock (GPa). Km (float): Bulk modulus of the matrix (GPa). phi (float): Porosity of the rock. Sat (float): Saturation level of the fluid in the pores.

Returns: float: Saturated bulk modulus (GPa).

PyHydroGeophysX.petrophysics.velocity_models.velDEM(phi, Km, Gm, rho_b, Sat, alpha)[source]

Calculate effective bulk modulus (Keff), shear modulus (Geff), and P-wave velocity (Vp) for a rock with varying porosity (phi) based on the DEM model, taking into account the saturation (Sat) and the crack aspect ratio (alpha).

Parameters: phi (np.array): Array of porosities. Km (float): Initial bulk modulus of the material (GPa). Gm (float): Initial shear modulus of the material (GPa). rho_b (float): Density of the solid phase (kg/m^3). Sat (float): Saturation level of the fluid in the cracks (0 to 1, where 1 is fully saturated). alpha (float): Crack aspect ratio.

Returns: Keff1 (np.array): Effective bulk modulus for each porosity value (GPa). Geff1 (np.array): Effective shear modulus for each porosity value (GPa). Vp (np.array): P-wave velocity for each porosity value (m/s).

PyHydroGeophysX.petrophysics.velocity_models.vel_porous(phi, Km, Gm, rho_b, Sat, depth=1)[source]

Calculate P-wave velocity (Vp) for a rock with varying porosity (phi) based on the Hertz-Mindlin model and Hashin-Shtrikman bounds, taking into account the saturation (Sat).

Parameters: phi (np.array): Array of porosities. Km (float): Bulk modulus of the solid phase (GPa). Gm (float): Shear modulus of the solid phase (GPa). rho_b (float): Density of the solid phase (kg/m^3). Sat (float): Saturation level of the fluid in the pores (0 to 1, where 1 is fully saturated). depth (float): depth for pressure estimation (m)

Returns: Vp_h (np.array): P-wave velocity for each porosity value (upper Hashin-Shtrikman bound) (m/s). Vp_l (np.array): P-wave velocity for each porosity value (lower Hashin-Shtrikman bound) (m/s).

Module contents