PyFVTool: Python toolbox for the finite-volume method
This is a Python implementation of A. A. Eftekhari's Matlab/Octave FVM solver FVTool heavily inspired by FiPy albeit with a fraction of FiPy features. The boundary conditions, however, are much easier (and perhaps more consistent) to implement in PyFVTool.
This package can discretize and solve the conservative form of transient convection-diffusion-reaction equation(s) with variable velocity field/diffusion coefficients:
with the following general form of boundary conditions (by specifying constants a
, b
, and c
):
The finite-volume discretization schemes include:
- 1D, 2D and 3D Cartesian and Cylindrical grids (only 1D spherical)
- Second order (central difference) diffusion terms
- Second order (central difference), first order (upwind), and total variation diminishing (TVD) for advection terms
- Constant and linear source terms
- Backward and forward Euler for transient terms
- Dirichlet, Neumann, Robin, and periodic boundary conditions
- (Relatively) easy linearization of nonlinear PDEs
- Averaging methods (linear, arithmetic, geometric, harmonic, upwind, TVD)
- Divergence and gradient terms
The code is under active development. Preliminary simulations match analytical solutions. More validation is under way, and the use of this PyFVTool toolbox in ongoing research projects will further consolidate the code and verify its validity.
There are simple tricks to use the code for systems of nonlinear partial differential equations that take the form of advection-diffusion-reaction equations. There is not much documentation for the code yet (help wanted!) but the example folder is the best place to start. If you know the topic, there is a good chance you will never need any documentations.
Installation
An important feature of PyFVTool is that it is 'pure scientific Python' (i.e. it needs only Python and the standard scientific computing libraries numpy
, scipy
and matplotlib
to run). Further optional dependencies may appear in the future, e.g., for increasing the computational speed via optimised numerical libraries, but these will remain optional.
For now, install it directly from the repo using pip. You will need Python 3.8
or higher and numpy
, scipy
, and matplotlib
:
pip install git+https://github.com/simulkade/PyFVTool.git
Soon, it will be possible to install from conda-forge and Python Package Index too.
Example
Here is a simple example of 1D transient diffusion equation:
from pyfvtool import *
# Solving a 1D diffusion equation with a fixed concentration
# at the left boundary and a closed boundary on the right side
Nx = 20 # number of finite volume cells
Lx = 1.0 # [m] length of the domain
c_left = 1.0 # left boundary concentration
c_init = 0.0 # initial concentration
D_val = 1e-5 # diffusion coefficient (gas phase)
t_simulation = 3600.0 # [s] simulation time
dt = 60.0 # [s] time step
m1 = createMesh1D(Nx, Lx) # mesh object
bc = createBC(m1) # Neumann boundary condition by default
# switch the left boundary to Dirichlet: fixed concentration
bc.left.a[:] = 0.0
bc.left.b[:] = 1.0
bc.left.c[:] = c_left
# create a cell variable with initial concentration
c_old = createCellVariable(m1, c_init, bc)
# assign diffusivity to cells
D_cell = createCellVariable(m1, D_val)
D_face = geometricMean(D_cell) # average value of diffusivity at the interfaces between cells
# Discretization
Mdiff = diffusionTerm(D_face)
Mbc, RHSbc = boundaryConditionTerm(bc)
# time loop
t = 0
while t<t_simulation:
t+=dt
Mt, RHSt = transientTerm(c_old, dt, 1.0)
c_new = solvePDE(m1, Mt-Mdiff+Mbc, RHSbc+RHSt)
c_old.update_value(c_new)
visualizeCells(c_old)