Berea Sandstone Simulation Using PoreSpy and OpenPNM

The example explains effective permeabilty calculations using PoreSpy and OpenPNM software. The simulation is performed on X-ray tomography image of BereaSandstone. The calculated effective permeablity value can compared with value report in Dong et al.

Start by importing the necessary packages


In [1]:
import os
import imageio
import scipy as sp
import numpy as np
import openpnm as op
import porespy as ps
import matplotlib.pyplot as plt
np.set_printoptions(precision=4)
np.random.seed(10)
%matplotlib inline

Load BreaSandstone Image file

Give path to image file and load the image. Please note image should be binarized or in boolean format before performing next steps.


In [2]:
path = '../fixtures/ICL-Sandstone(Berea)/'
file_format = '.tif'
file_name = 'Berea'
file = file_name + file_format
fetch_file = os.path.join(path, file)
im = imageio.mimread(fetch_file)
im = ~np.array(im, dtype=bool)

Confirm image and check image porosity

Be patient, this might take ~30 seconds (depending on your CPU)


In [3]:
#NBVAL_IGNORE_OUTPUT
fig, ax = plt.subplots(1, 3, figsize=(12,5))
ax[0].imshow(im[:, :, 100]);
ax[1].imshow(ps.visualization.show_3D(im));
ax[2].imshow(ps.visualization.sem(im));
ax[0].set_title("Slice No. 100 View");
ax[1].set_title("3D Sketch");
ax[2].set_title("SEM View");



In [4]:
print(ps.metrics.porosity(im))


0.19645303125

Extract pore network using SNOW algorithm in PoreSpy

The SNOW algorithm (an accronym for Sub-Network from an Oversegmented Watershed) was presented by Gostick. The algorithm was used to extract pore network from BereaSandstone image.


In [5]:
#NBVAL_IGNORE_OUTPUT
resolution = 5.345e-6
net = ps.networks.snow(im=im, voxel_size=resolution)


____________________________________________________________
Beginning SNOW Algorithm
Converting supplied image (im) to boolean
Peforming Distance Transform
Applying Gaussian blur with sigma = 0.4
Initial number of peaks:  10878
Peaks after trimming saddle points:  4315
Peaks after trimming nearby peaks:  4315
____________________________________________________________
Extracting pore and throat information from image
100%|██████████| 5478/5478 [00:50<00:00, 109.51it/s]

Import network in OpenPNM

The output from the SNOW algorithm above is a plain python dictionary containing all the extracted pore-scale data, but it is NOT yet an OpenPNM network. We need to create an empty network in OpenPNM, then populate it with the data from SNOW:


In [6]:
ws = op.Workspace()
proj = op.Project()
pn = op.network.GenericNetwork(name=file_name, project=proj)
pn.update(net)  # Fills 'pn' with data from 'net'

Now we can print the network to see how the transferred worked:


In [7]:
print(pn)


――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――
openpnm.network.GenericNetwork : Berea
――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――
#     Properties                                    Valid Values
――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――
1     pore.area                                      5478 / 5478 
2     pore.centroid                                  5478 / 5478 
3     pore.coords                                    5478 / 5478 
4     pore.diameter                                  5478 / 5478 
5     pore.equivalent_diameter                       5478 / 5478 
6     pore.extended_diameter                         5478 / 5478 
7     pore.inscribed_diameter                        5478 / 5478 
8     pore.label                                     5478 / 5478 
9     pore.surface_area                              5478 / 5478 
10    pore.volume                                    5478 / 5478 
11    throat.area                                    9523 / 9523 
12    throat.centroid                                9523 / 9523 
13    throat.conduit_lengths.pore1                   9523 / 9523 
14    throat.conduit_lengths.pore2                   9523 / 9523 
15    throat.conduit_lengths.throat                  9523 / 9523 
16    throat.conns                                   9523 / 9523 
17    throat.diameter                                9523 / 9523 
18    throat.direct_length                           9523 / 9523 
19    throat.endpoints.head                          9523 / 9523 
20    throat.endpoints.tail                          9523 / 9523 
21    throat.equivalent_diameter                     9523 / 9523 
22    throat.inscribed_diameter                      9523 / 9523 
23    throat.length                                  9523 / 9523 
24    throat.perimeter                               9523 / 9523 
25    throat.total_length                            9523 / 9523 
26    throat.volume                                  9523 / 9523 
――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――
#     Labels                                        Assigned Locations
――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――
1     pore.all                                      5478      
2     pore.back                                     187       
3     pore.bottom                                   196       
4     pore.boundary                                 1163      
5     pore.front                                    189       
6     pore.internal                                 4315      
7     pore.left                                     179       
8     pore.right                                    206       
9     pore.top                                      206       
10    throat.all                                    9523      
11    throat.boundary                               1160      
12    throat.internal                               8363      
――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――

Check network health

Remove isolated pores or cluster of pores from the network by checking it network health. Make sure ALL keys in network health functions have no value.


In [8]:
h = pn.check_network_health()
op.topotools.trim(network=pn, pores=h['trim_pores'])
h = pn.check_network_health()
print(h)


――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――
key                                 value
――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――
disconnected_clusters               []
isolated_pores                      []
trim_pores                          []
duplicate_throats                   []
bidirectional_throats               []
headless_throats                    []
looped_throats                      []
――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――――

Assign pore scale model

Assign conduit shape to calculate flow from pore i to pore j. In this simulation spherical pore and cylinderical throats geometry is considered.


In [9]:
mod = op.models.geometry.throat_endpoints.spherical_pores
pn.add_model(propname='throat.endpoints', model=mod)
mod = op.models.geometry.throat_length.conduit_lengths
pn.add_model(propname='throat.conduit_lengths', model=mod)
mod = op.models.geometry.pore_area.sphere
pn.add_model(propname='pore.area', model=mod)

Assign geometry


In [10]:
geo = op.geometry.GenericGeometry(network=pn, pores=pn.Ps, throats=pn.Ts)

Assign phase

In this example air is considered as fluid passing through porous channels.


In [11]:
air = op.phases.Air(network=pn)

Assign physics


In [12]:
phys_air = op.physics.Standard(network=pn, phase=air, geometry=geo)

Assign Algorithm and boundary conditions

Select stokes flow algorithm for simulation and assign dirichlet boundary conditions in top and bottom faces of the network.


In [13]:
perm = op.algorithms.StokesFlow(network=pn, project=proj)
perm.setup(phase=air)
perm.set_value_BC(pores=pn.pores('top'), values=0)
perm.set_value_BC(pores=pn.pores('bottom'), values=101325)
perm.run()
air.update(perm.results())

Calculate effective permeability

Caclulate effective permeablity using hagen poiseuille equation. Use cross section area and flow length manually from image dimension.


In [14]:
resolution = 5.345e-6
Q = perm.rate(pores=pn.pores('bottom'), mode='group')
A = (im.shape[0] * im.shape[1]) * resolution**2
L = im.shape[2] * resolution
mu = air['pore.viscosity'].max()
delta_P = 101325 - 0
K = Q * L * mu / (A * delta_P)
print('The value of K is:', K/0.98e-12*1000, 'mD')


The value of K is: [944.4486] mD