API reference

This section documents the public API of PyFVTool and is generated automatically from the source code docstrings.

To do: API documentation

  • improve docstrings

  • re-verify the automodule directives (in docs/api/index.md)

  • docs/api/index.md is just a plain Markdown file, it can be edited freely. Add the missing modules and restructure however makes sense.

Meshes

class pyfvtool.mesh.CellLocation(_x: ndarray, _y: ndarray, _z: ndarray, coordlabels: dict)[source]

Bases: CellProp

class pyfvtool.mesh.CellProp(_x: ndarray, _y: ndarray, _z: ndarray, coordlabels: dict)[source]

Bases: object

property phi
property r
property theta
property x
property y
property z
class pyfvtool.mesh.CellSize(_x: ndarray, _y: ndarray, _z: ndarray, coordlabels: dict)[source]

Bases: CellProp

class pyfvtool.mesh.CylindricalGrid1D(Nr: int, Lr: float)[source]
class pyfvtool.mesh.CylindricalGrid1D(face_locationsR: ndarray)
class pyfvtool.mesh.CylindricalGrid1D(dims, cellsize, cellcenters, facecenters, corners, edges)

Bases: Grid1D

Mesh based on a 1D cylindrical grid (r)

This class can be instantiated in different ways: from a list of cell face locations or from the number of cells and domain length.

Instantiation Options:

  • CylindricalGrid1D(Nr, Lr)

  • CylindricalGrid1D(face_locationsR)

param CylindricalGrid1D(Nr:
Nrint

Number of cells in the r direction.

Lrfloat

Length of the domain in the r direction.

param Lr):
Nrint

Number of cells in the r direction.

Lrfloat

Length of the domain in the r direction.

param CylindricalGrid1D(face_locationsR):
face_locationsRndarray

Locations of the cell faces in the r direction.

Examples

>>> import numpy as np
>>> from pyfvtool import CylindricalGrid1D
>>> mesh = CylindricalGrid1D(10, 10.0)
>>> print(mesh)
class pyfvtool.mesh.CylindricalGrid2D(Nr: int, Nz: int, Lr: float, Lz: float)[source]
class pyfvtool.mesh.CylindricalGrid2D(face_locationsR: ndarray, face_locationsZ: ndarray)
class pyfvtool.mesh.CylindricalGrid2D(dims, cellsize, cellcenters, facecenters, corners, edges)

Bases: Grid2D

Mesh based on a 2D cylindrical grid (r, z)

This class can be instantiated in different ways: from a list of cell face locations or from the number of cells and domain length.

Instantiation Options:

  • CylindricalGrid2D(Nr, Nz, Lr, Lz)

  • CylindricalGrid2D(face_locationsR, face_locationsZ)

param CylindricalGrid2D(Nr:
Nrint

Number of cells in the r direction.

Nzint

Number of cells in the z direction.

Lrfloat

Length of the domain in the r direction.

Lzfloat

Length of the domain in the z direction.

param Nz:
Nrint

Number of cells in the r direction.

Nzint

Number of cells in the z direction.

Lrfloat

Length of the domain in the r direction.

Lzfloat

Length of the domain in the z direction.

param Lr:
Nrint

Number of cells in the r direction.

Nzint

Number of cells in the z direction.

Lrfloat

Length of the domain in the r direction.

Lzfloat

Length of the domain in the z direction.

param Lz):
Nrint

Number of cells in the r direction.

Nzint

Number of cells in the z direction.

Lrfloat

Length of the domain in the r direction.

Lzfloat

Length of the domain in the z direction.

param CylindricalGrid2D(face_locationsR:
face_locationsRndarray

Locations of the cell faces in the r direction.

face_locationsZndarray

Locations of the cell faces in the z direction.

param face_locationsZ):
face_locationsRndarray

Locations of the cell faces in the r direction.

face_locationsZndarray

Locations of the cell faces in the z direction.

returns:

A 2D cylindrical mesh object.

rtype:

CylindricalGrid2D

Examples

>>> import numpy as np
>>> from pyfvtool import CylindricalGrid2D
>>> mesh = CylindricalGrid2D(10, 10, 10.0, 10.0)
>>> print(mesh)
class pyfvtool.mesh.CylindricalGrid3D(Nr: int, Ntheta: int, Nz: int, Lr: float, Ltheta: float, Lz: float)[source]
class pyfvtool.mesh.CylindricalGrid3D(face_locationsR: ndarray, face_locationsTheta: ndarray, face_locationsZ: ndarray)
class pyfvtool.mesh.CylindricalGrid3D(dims, cellsize, cellcenters, facecenters, corners, edges)

Bases: Grid3D

Mesh based on a 3D cylindrical grid (r, theta, z)

This class can be instantiated in different ways: from a list of cell face locations or from the number of cells and domain length.

Instantiation Options:

  • CylindricalGrid3D(Nr, Ntheta, Nz, Lr, Ltheta, Lz)

  • CylindricalGrid3D(face_locationsR, face_locationsTheta, face_locationsZ)

param CylindricalGrid3D(Nr:
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nzint

Number of cells in the z direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lzfloat

Length of the domain in the z direction.

param Ntheta:
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nzint

Number of cells in the z direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lzfloat

Length of the domain in the z direction.

param Nz:
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nzint

Number of cells in the z direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lzfloat

Length of the domain in the z direction.

param Lr:
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nzint

Number of cells in the z direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lzfloat

Length of the domain in the z direction.

param Ltheta:
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nzint

Number of cells in the z direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lzfloat

Length of the domain in the z direction.

param Lz):
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nzint

Number of cells in the z direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lzfloat

Length of the domain in the z direction.

CylindricalGrid3D(face_locationsR, face_locationsTheta, face_locationsZ)
face_locationsRndarray

Locations of the cell faces in the r direction.

face_locationsThetandarray

Locations of the cell faces in the theta direction.

face_locationsZndarray

Locations of the cell faces in the z direction.

Examples

>>> import numpy as np
>>> from pyfvtool import CylindricalGrid3D
>>> mesh = CylindricalGrid3D(10, 10, 10, 10.0, 10.0, 10.0)
>>> print(mesh)
class pyfvtool.mesh.FaceLocation(_x: ndarray, _y: ndarray, _z: ndarray, coordlabels: dict)[source]

Bases: CellProp

class pyfvtool.mesh.Grid1D(Nx: int, Lx: float)[source]
class pyfvtool.mesh.Grid1D(face_locations: ndarray)
class pyfvtool.mesh.Grid1D(dims, cellsize, cellcenters, facecenters, corners, edges)

Bases: MeshStructure

Mesh based on a 1D Cartesian grid (x)

This class can be instantiated in different ways: from a list of cell face locations or from the number of cells and domain length.

Instantiation Options:

  • Grid1D(Nx, Lx)

  • Grid1D(face_locationsX)

param Grid1D(Nx:
Nxint

Number of cells in the x direction.

Lxfloat

Length of the domain in the x direction.

param Lx):
Nxint

Number of cells in the x direction.

Lxfloat

Length of the domain in the x direction.

param Grid1D(face_locationsX):
face_locationsXndarray

Locations of the cell faces in the x direction.

Examples

>>> import numpy as np
>>> from pyfvtool import Grid1D
>>> mesh = Grid1D(10, 10.0)
>>> print(mesh)
cell_numbers()[source]
class pyfvtool.mesh.Grid2D(Nx: int, Ny: int, Lx: float, Ly: float)[source]
class pyfvtool.mesh.Grid2D(face_locationsX: ndarray, face_locationsY: ndarray)
class pyfvtool.mesh.Grid2D(dims, cellsize, cellcenters, facecenters, corners, edges)

Bases: MeshStructure

Mesh based on a 2D Cartesian grid (x, y)

This class can be instantiated in different ways: from a list of cell face locations or from the number of cells and domain length.

Instantiation Options:

  • Grid2D(Nx, Ny, Lx, Ly)

  • Grid2D(face_locationsX, face_locationsY)

param Grid2D(Nx:
Nxint

Number of cells in the x direction.

Nyint

Number of cells in the y direction.

Lxfloat

Length of the domain in the x direction.

Lyfloat

Length of the domain in the y direction.

param Ny:
Nxint

Number of cells in the x direction.

Nyint

Number of cells in the y direction.

Lxfloat

Length of the domain in the x direction.

Lyfloat

Length of the domain in the y direction.

param Lx:
Nxint

Number of cells in the x direction.

Nyint

Number of cells in the y direction.

Lxfloat

Length of the domain in the x direction.

Lyfloat

Length of the domain in the y direction.

param Ly):
Nxint

Number of cells in the x direction.

Nyint

Number of cells in the y direction.

Lxfloat

Length of the domain in the x direction.

Lyfloat

Length of the domain in the y direction.

param Grid2D(face_locationsX:
face_locationsXndarray

Locations of the cell faces in the x direction.

face_locationsYndarray

Locations of the cell faces in the y direction.

param face_locationsY):
face_locationsXndarray

Locations of the cell faces in the x direction.

face_locationsYndarray

Locations of the cell faces in the y direction.

Examples

>>> import numpy as np
>>> from pyfvtool import Grid2D
>>> mesh = Grid2D(10, 10, 10.0, 10.0)
>>> print(mesh)
cell_numbers()[source]
class pyfvtool.mesh.Grid3D(Nx: int, Ny: int, Nz: int, Lx: float, Ly: float, Lz: float)[source]
class pyfvtool.mesh.Grid3D(dims, cellsize, cellcenters, facecenters, corners, edges)

Bases: MeshStructure

Mesh based on a 3D Cartesian grid (x, y, z)

This class can be instantiated in different ways: from a list of cell face locations or from the number of cells and domain length. There are multiple overloaded ‘__init__’ methods available to provide flexibility in instantiation.

Instantiation Options:

  • Grid3D(Nx, Ny, Nz, Lx, Ly, Lz)

  • Grid3D(face_locationsX, face_locationsY, face_locationsZ)

param Grid3D(Nx:
Nxint

Number of cells in the x direction.

Nyint

Number of cells in the y direction.

Nzint

Number of cells in the z direction.

Lxfloat

Length of the domain in the x direction.

Lyfloat

Length of the domain in the y direction.

Lzfloat

Length of the domain in the z direction.

param Ny:
Nxint

Number of cells in the x direction.

Nyint

Number of cells in the y direction.

Nzint

Number of cells in the z direction.

Lxfloat

Length of the domain in the x direction.

Lyfloat

Length of the domain in the y direction.

Lzfloat

Length of the domain in the z direction.

param Nz:
Nxint

Number of cells in the x direction.

Nyint

Number of cells in the y direction.

Nzint

Number of cells in the z direction.

Lxfloat

Length of the domain in the x direction.

Lyfloat

Length of the domain in the y direction.

Lzfloat

Length of the domain in the z direction.

param Lx:
Nxint

Number of cells in the x direction.

Nyint

Number of cells in the y direction.

Nzint

Number of cells in the z direction.

Lxfloat

Length of the domain in the x direction.

Lyfloat

Length of the domain in the y direction.

Lzfloat

Length of the domain in the z direction.

param Ly:
Nxint

Number of cells in the x direction.

Nyint

Number of cells in the y direction.

Nzint

Number of cells in the z direction.

Lxfloat

Length of the domain in the x direction.

Lyfloat

Length of the domain in the y direction.

Lzfloat

Length of the domain in the z direction.

param Lz):
Nxint

Number of cells in the x direction.

Nyint

Number of cells in the y direction.

Nzint

Number of cells in the z direction.

Lxfloat

Length of the domain in the x direction.

Lyfloat

Length of the domain in the y direction.

Lzfloat

Length of the domain in the z direction.

param Grid3D(face_locationsX:
face_locationsXndarray

Locations of the cell faces in the x direction.

face_locationsYndarray

Locations of the cell faces in the y direction.

face_locationsZndarray

Locations of the cell faces in the z direction.

param face_locationsY:
face_locationsXndarray

Locations of the cell faces in the x direction.

face_locationsYndarray

Locations of the cell faces in the y direction.

face_locationsZndarray

Locations of the cell faces in the z direction.

param face_locationsZ):
face_locationsXndarray

Locations of the cell faces in the x direction.

face_locationsYndarray

Locations of the cell faces in the y direction.

face_locationsZndarray

Locations of the cell faces in the z direction.

Examples

>>> import numpy as np
>>> from pyfvtool import Grid3D
>>> mesh = Grid3D(10, 10, 10, 10.0, 10.0, 10.0)
>>> print(mesh)
cell_numbers()[source]
class pyfvtool.mesh.MeshStructure(dims, cellsize, cellcenters, facecenters, corners, edges)[source]

Bases: object

property cellvolume
class pyfvtool.mesh.PolarGrid2D(Nr: int, Ntheta: int, Lr: float, Ltheta: float)[source]
class pyfvtool.mesh.PolarGrid2D(face_locationsR: ndarray, face_locationsTheta: ndarray)
class pyfvtool.mesh.PolarGrid2D(dims, cellsize, cellcenters, facecenters, corners, edges)

Bases: Grid2D

Mesh based on a 2D polar grid (r, theta)

This class can be instantiated in different ways: from a list of cell face locations or from the number of cells and domain length.

Instantiation Options:

  • PolarGrid2D(Nr, Ntheta, Lr, Ltheta)

  • PolarGrid2D(face_locationsR, face_locationsTheta)

Parameters:

PolarGrid2D(Nr, Ntheta, Lr, Ltheta):
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

PolarGrid2D(face_locationsR, face_locationsTheta):
face_locationsRndarray

Locations of the cell faces in the r direction.

face_locationsThetandarray

Locations of the cell faces in the theta direction.

Examples:

>>> import numpy as np
>>> from pyfvtool import PolarGrid2D
>>> mesh = PolarGrid2D(10, 10, 10.0, 10.0)
>>> print(mesh)
class pyfvtool.mesh.SphericalGrid1D(Nr: int, Lr: float)[source]
class pyfvtool.mesh.SphericalGrid1D(face_locationsR: ndarray)
class pyfvtool.mesh.SphericalGrid1D(dims, cellsize, cellcenters, facecenters, corners, edges)

Bases: Grid1D

Mesh based on a 1D spherical grid (r)

This class can be instantiated in different ways: from a list of cell face locations or from the number of cells and domain length.

Instantiation Options:

  • SphericalGrid1D(Nr, Lr)

  • SphericalGrid1D(face_locationsR)

param SphericalGrid1D(Nr:
Nrint

Number of cells in the r direction.

Lrfloat

Length of the domain in the r direction.

param Lr):
Nrint

Number of cells in the r direction.

Lrfloat

Length of the domain in the r direction.

param SphericalGrid1D(face_locationsR):
face_locationsRndarray

Locations of the cell faces in the r direction.

Examples

>>> import numpy as np
>>> from pyfvtool import SphericalGrid1D
>>> mesh = SphericalGrid1D(10, 10.0)
>>> print(mesh)
class pyfvtool.mesh.SphericalGrid3D(Nr: int, Ntheta: int, Nphi: int, Lr: float, Ltheta: float, Lphi: float)[source]
class pyfvtool.mesh.SphericalGrid3D(face_locationsR: ndarray, face_locationsTheta: ndarray, face_locationsPhi: ndarray)
class pyfvtool.mesh.SphericalGrid3D(dims, cellsize, cellcenters, facecenters, corners, edges)

Bases: Grid3D

Mesh based on a 3D spherical grid (r, theta, phi)

Create a SphericalGrid3D object from a list of cell face locations or from the number of cells and domain length.

Instantiation Options:

  • SphericalGrid3D(Nr, Ntheta, Nphi, Lr, Ltheta, Lphi)

  • SphericalGrid3D(face_locationsR, face_locationsTheta, face_locationsPhi)

param SphericalGrid3D(Nr:
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nphiint

Number of cells in the phi direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lphifloat

Length of the domain in the phi direction.

param Ntheta:
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nphiint

Number of cells in the phi direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lphifloat

Length of the domain in the phi direction.

param Nphi:
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nphiint

Number of cells in the phi direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lphifloat

Length of the domain in the phi direction.

param Lr:
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nphiint

Number of cells in the phi direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lphifloat

Length of the domain in the phi direction.

param Ltheta:
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nphiint

Number of cells in the phi direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lphifloat

Length of the domain in the phi direction.

param Lphi):
Nrint

Number of cells in the r direction.

Nthetaint

Number of cells in the theta direction.

Nphiint

Number of cells in the phi direction.

Lrfloat

Length of the domain in the r direction.

Lthetafloat

Length of the domain in the theta direction.

Lphifloat

Length of the domain in the phi direction.

SphericalGrid3D(face_locationsR, face_locationsTheta, face_locationsPhi)
face_locationsRndarray

Locations of the cell faces in the r direction.

face_locationsThetandarray

Locations of the cell faces in the theta direction.

face_locationsPhindarray

Locations of the cell faces in the phi direction.

Examples

>>> import numpy as np
>>> from pyfvtool import SphericalGrid3D
>>> mesh = SphericalGrid3D(10, 10, 10, 10.0, 10.0, 10.0)
>>> print(mesh)

Cell and face variables

Cell variables

class pyfvtool.cell.CellVariable(mesh_struct: MeshStructure, cell_value: ndarray, BC: BoundaryConditionsBase)[source]
class pyfvtool.cell.CellVariable(mesh_struct: MeshStructure, cell_value: ndarray)
class pyfvtool.cell.CellVariable(mesh_struct: MeshStructure, cell_value: float, BC: BoundaryConditionsBase)
class pyfvtool.cell.CellVariable(mesh_struct: MeshStructure, cell_value: float)

Bases: object

apply_BCs()[source]

Update ghost cells according to the boundary conditions and the internal (inner) cell values.

In standard use, i.e. when solving equations with solvePDE(), it is NOT necessary to call this method, as the updates are handled automatically.

In certain special ‘expert-level’ cases, specific calls to this method may be necessary. For instance, when the CellVariable is used before calling solvePDE(). Or when when working with the function solveMatrixPDE() instead of solvePDE().

In general, superfluous calls to apply_BCs() will not hurt.

Returns:

  • None.

  • The modified attribute of the CellVariableBCs is reset, as well as

  • the modified attribute of the CellVariable.value

property cellcenters
property cellvolume
copy()[source]

Create a copy of the CellVariable

Returns:

Copy of the CellVariable.

Return type:

CellVariable

domainIntegral() float[source]

Calculate the finite-volume integral of a CellVariable over entire domain

The finite-volume integral over the entire mesh domain gives the total amount of CellVariable present in the system. Calculation of this integral is useful for checking conservation of the quantity concerned (in case of ‘no-flux’ BCs), or for monitoring its evolution due to exchanges via the boundaries (other BCs) or to the presence of source terms.

May later become a built-in method of the CellVariable class, but for now this implementation as a function is chosen for consistency with FVTool.

Returns:

Total finite-volume integral over entire domain.

Return type:

float

plotprofile()[source]

Create a ‘profile’ of a cell variable for plotting, export, etc.

‘Plot profiles’ are sets of arrays that contain the axes’ coordinates, cell values, including the values at the boundaries.

For 2D and 3D visualization, it is perhaps better to use only the actual cell values (e.g. for a false color map à la plt.pcolormesh), and not include the values at the boundaries.

1D meshes

This generates a pair of vectors containing the abscissa and ordinates for plotting the values of the cell variable over the entire calculation domain. It includes the values at the outer faces of the domain, by taking into account the values of the ghost cells.

returns:
  • x (np.ndarray) – x (or r) coordinates.

  • phi0 (np.ndarray) – Value of the CellVariables at those points.

2D meshes

This generates a set of vectors containing the (x, y) or (r, z) coordinates and the values of the cell variable at those coordinates for plotting the values of the cell variable over the entire calculation domain. It includes the values at the outer faces of the domain, by taking into account the values of the ghost cells.

returns:
  • x (np.ndarray) – x (or r) coordinates.

  • y (np.ndarray) – y (or z) coordinates.

  • phi0 (np.ndarray) – Value of the CellVariables at those coordinates.

3D meshes

This generates a set of vectors containing the (x, y, z) or (r, theta, z) coordinates and the values of the cell variable at those coordinates for plotting the values of the cell variable over the entire calculation domain. It includes the values at the outer faces of the domain, by taking into account the values of the ghost cells.

returns:
  • x (np.ndarray) – x (or r) coordinates.

  • y (np.ndarray) – y (or theta) coordinates.

  • z (np.ndarray) – z coordinates.

  • phi0 (np.ndarray) – Value of the CellVariables at those coordinates.

update_value(new_cell)[source]
property value
pyfvtool.cell.cellLocations(m: MeshStructure)[source]

this function returns the location of the cell centers as cell variables.

It can later be used in defining properties that are variable in space.

Incompletely tested, and there may be other, more direct, ways to calculate properties that vary in space, e.g. using cellcenters directly?

Parameters:

m ({MeshStructure object}) – Domain of the problem

Returns:

  • X ({CellVariable object}) – Node x-positions

  • Y ({CellVariable object}) – Node y-positions

  • Z ({CellVariable object}) – Node z-positions

See also

faceLocations

Notes

Examples

>>>
pyfvtool.cell.celleval(f, *args)[source]
pyfvtool.cell.funceval(f, *args)[source]

Face variables

class pyfvtool.face.FaceVariable(mesh: MeshStructure, faceval: float)[source]
class pyfvtool.face.FaceVariable(mesh: MeshStructure, faceval: ndarray)
class pyfvtool.face.FaceVariable(mesh: MeshStructure, _xvalue: ndarray, _yvalue: ndarray, _zvalue: ndarray)

Bases: object

Face variable class

Create a FaceVariable for the given mesh with the given value. Examples: >>> import pyfvtool as pf >>> m = pf.Grid1D(10, 1.0) >>> f = pf.FaceVariable(m, 1.0)

property phivalue
property rvalue
property thetavalue
property xvalue
property yvalue
property zvalue
pyfvtool.face.faceLocations(m: MeshStructure)[source]

this function returns the location of the cell faces as face variables.

It can later be used for the calculation of face variables as a function of location

Incompletely tested

Parameters:

m ({MeshStructure object}) – Domain of the problem

Returns:

  • X ({FaceVariable object}) – Edge x-positions

  • Y ({FaceVariable object}) – Edge y-positions

  • Z ({FaceVariable object}) – Edge z-positions

See also

cellLocations

Notes

Examples

>>>
pyfvtool.face.faceeval(f, *args)[source]

Evaluate a function f on a FaceVariable. Examples: >>> import pyfvtool as pf >>> m = pf.Grid1D(10, 1.0) >>> f = pf.FaceVariable(m, 1.0) >>> g = pf.faceeval(lambda x: x**2, f)

Boundary conditions

pyfvtool.boundary.BoundaryConditions(mesh: MeshStructure)[source]

Create an object for holding boundary conditions (factory function)

Parameters:

mesh (MeshStructure) – Mesh used for space discretization.

Returns:

Object holding all information defining boundary conditions.

Return type:

Subclass instance of BoundaryConditionsBase class

class pyfvtool.boundary.BoundaryConditions1D(mesh: Grid1D)[source]

Bases: BoundaryConditionsBase

class pyfvtool.boundary.BoundaryConditions2D(mesh: Grid2D)[source]

Bases: BoundaryConditionsBase

class pyfvtool.boundary.BoundaryConditions3D(mesh: Grid3D)[source]

Bases: BoundaryConditionsBase

class pyfvtool.boundary.BoundaryConditionsBase(mesh: MeshStructure, left: BoundaryFace, right: BoundaryFace, bottom: BoundaryFace, top: BoundaryFace, back: BoundaryFace, front: BoundaryFace)[source]

Bases: object

Base class for the boundary conditions

mesh

mesh structure

Type:

MeshStructure

left

boundary condition for the left face

Type:

BoundaryFace

right

boundary condition for the right face

Type:

BoundaryFace

bottom

boundary condition for the bottom face

Type:

BoundaryFace

top

boundary condition for the top face

Type:

BoundaryFace

back

boundary condition for the back face

Type:

BoundaryFace

front

boundary condition for the back face

Type:

BoundaryFace

property modified

True if any of the BoundaryFace conditions has changed since last apply_BCs()

class pyfvtool.boundary.BoundaryFace(a: ndarray, b: ndarray, c: ndarray, periodic=False)[source]

Bases: object

Class describing the boundary condition of a single face

property a

Coefficient array a of the boundary condition.

property b

Coefficient array b of the boundary condition.

property c

Coefficient array c of the boundary condition.

property modified

Boolean to indicate if boundary conditions were modified since last update.

property periodic

Boolean. If True, use periodic boundary condition.

pyfvtool.boundary.boundaryConditionsTerm(BC)[source]

Generate the terms of the matrix equation representing the boundary conditions

Parameters:

BC (instance of subclass of BoundaryConditions) – Boundary conditions

Returns:

  • BCMatrix (csr_array) – Matrix elements representing the boundary conditions

  • BCRHS (numpy.ndarray (1D)) – RHS column vector representing the boundary conditions

pyfvtool.boundary.boundaryConditionsTerm1D(BC: BoundaryConditions1D)[source]
pyfvtool.boundary.boundaryConditionsTerm2D(BC: BoundaryConditions2D)[source]
pyfvtool.boundary.boundaryConditionsTerm3D(BC: BoundaryConditions3D)[source]
pyfvtool.boundary.boundaryConditionsTermCylindrical3D(BC: BoundaryConditions3D)[source]
pyfvtool.boundary.boundaryConditionsTermPolar2D(BC: BoundaryConditions2D)[source]
pyfvtool.boundary.boundaryConditionsTermSpherical3D(BC: BoundaryConditions3D)[source]
pyfvtool.boundary.cellValuesWithBoundaries(phi, BC) ndarray[source]

Return all cell values with the boundary values calculated given the boundary conditions.

Parameters:
  • phi (array_like) – internal cell values

  • BC (BoundaryCondition) – Boundary condition object

Returns:

phiBC – All cell values (internal and boundary), including boundaries

Return type:

array_like

pyfvtool.boundary.cellValuesWithBoundaries1D(phi, BC)[source]

cellValuesWithBoundaries for 1D mesh

pyfvtool.boundary.cellValuesWithBoundaries2D(phi, BC)[source]

cellValuesWithBoundaries for 2D mesh

pyfvtool.boundary.cellValuesWithBoundaries3D(phi, BC)[source]

cellValuesWithBoundaries for 3D mesh

pyfvtool.boundary.cellValuesWithBoundariesCylindrical3D(phi, BC)[source]

cellValuesWithBoundaries for Cylindrical3D mesh

pyfvtool.boundary.cellValuesWithBoundariesPolar2D(phi, BC)[source]

cellValuesWithBoundaries for Polar2D mesh

pyfvtool.boundary.cellValuesWithBoundariesSpherical3D(phi, BC)[source]

cellValuesWithBoundaries for Spherical3D mesh

Discretization terms

Diffusion terms

pyfvtool.diffusion.diffusionTerm(D: FaceVariable) csr_array[source]

Builds the discretized diffusion term in matrix form for a given Cartesian mesh.

Parameters:

D (FaceVariable) – The diffusion coefficient.

Returns:

The discretized diffusion term in matrix form.

Return type:

csr_array

Examples

>>> import pyfvtool as pf
>>> m = pf.Grid1D(10, 1.0)
>>> D = pf.FaceVariable(m, 1.0)
>>> M = pf.diffusionTerm(D)
pyfvtool.diffusion.diffusionTerm1D(D: FaceVariable) csr_array[source]
pyfvtool.diffusion.diffusionTerm2D(D: FaceVariable) csr_array[source]
pyfvtool.diffusion.diffusionTerm3D(D: FaceVariable) csr_array[source]
pyfvtool.diffusion.diffusionTermCylindrical1D(D: FaceVariable) csr_array[source]
pyfvtool.diffusion.diffusionTermCylindrical2D(D: FaceVariable) csr_array[source]
pyfvtool.diffusion.diffusionTermCylindrical3D(D: FaceVariable) csr_array[source]
pyfvtool.diffusion.diffusionTermPolar2D(D: FaceVariable) csr_array[source]
pyfvtool.diffusion.diffusionTermSpherical1D(D: FaceVariable) csr_array[source]
pyfvtool.diffusion.diffusionTermSpherical3D(D: FaceVariable) csr_array[source]

Advection terms

pyfvtool.advection.convectionTVDupwindRHSTerm(u: FaceVariable, phi: CellVariable, FL, *args) ndarray[source]

Returns the TVD correction vector for the upwind convection term, \(\nabla \cdot (u \phi)\). The convection term is evaluated based on the TVD upwind scheme.

Parameters:

u (FaceVariable) – The velocity field.

Return type:

The discretized correction vector for the convection term.

Examples

>>> import pyfvtool as pf
>>> m = pf.Grid1D(10, 1.0)
>>> u = pf.FaceVariable(m, 1.0)
>>> phi = pf.CellVariable(m, 1.0)
>>> FL = pf.fluxLimiter('SUPERBEE')
>>> RHS = pf.convectionTVDupwindRHSTerm(u, phi, FL)
pyfvtool.advection.convectionTerm(u: FaceVariable) csr_array[source]

Returns the discretized convection term, \(\nabla \cdot (u \phi)\).

Parameters:

u (FaceVariable) – The velocity field.

Return type:

The discretized convection term.

Examples

>>> import pyfvtool as pf
>>> m = pf.Grid1D(10, 1.0)
>>> u = pf.FaceVariable(m, 1.0)
>>> M = pf.convectionTerm(u)
pyfvtool.advection.convectionTerm1D(u: FaceVariable)[source]
pyfvtool.advection.convectionTerm2D(u: FaceVariable)[source]
pyfvtool.advection.convectionTerm3D(u: FaceVariable)[source]
pyfvtool.advection.convectionTermCylindrical1D(u: FaceVariable)[source]
pyfvtool.advection.convectionTermCylindrical2D(u: FaceVariable)[source]
pyfvtool.advection.convectionTermCylindrical3D(u: FaceVariable)[source]
pyfvtool.advection.convectionTermPolar2D(u: FaceVariable)[source]
pyfvtool.advection.convectionTermSpherical1D(u: FaceVariable)[source]
pyfvtool.advection.convectionTermSpherical3D(u: FaceVariable)[source]
pyfvtool.advection.convectionTvdRHS1D(u: FaceVariable, phi: CellVariable, FL, *args)[source]
pyfvtool.advection.convectionTvdRHS2D(u: FaceVariable, phi: CellVariable, FL, *args)[source]
pyfvtool.advection.convectionTvdRHS3D(u: FaceVariable, phi: CellVariable, FL, *args)[source]
pyfvtool.advection.convectionTvdRHSCylindrical1D(u: FaceVariable, phi: CellVariable, FL, *args)[source]
pyfvtool.advection.convectionTvdRHSCylindrical2D(u: FaceVariable, phi: CellVariable, FL, *args)[source]
pyfvtool.advection.convectionTvdRHSCylindrical3D(u: FaceVariable, phi: CellVariable, FL, *args)[source]
pyfvtool.advection.convectionTvdRHSPolar2D(u: FaceVariable, phi: CellVariable, FL, *args)[source]
pyfvtool.advection.convectionTvdRHSSpherical1D(u: FaceVariable, phi: CellVariable, FL, *args)[source]
pyfvtool.advection.convectionTvdRHSSpherical3D(u: FaceVariable, phi: CellVariable, FL, *args)[source]
pyfvtool.advection.convectionUpwindTerm(u: FaceVariable, *args) csr_array[source]

Returns the discretized upwind convection term, \(\nabla \cdot (u \phi)\). The convection term is evaluated based on the upwind scheme.

Parameters:

u (FaceVariable) – The velocity field.

Return type:

The discretized upwind convection term.

Examples

>>> import pyfvtool as pf
>>> m = pf.Grid1D(10, 1.0)
>>> u = pf.FaceVariable(m, 1.0)
>>> M = pf.convectionUpwindTerm(u)
pyfvtool.advection.convectionUpwindTerm1D(u: FaceVariable, *args)[source]
pyfvtool.advection.convectionUpwindTerm2D(u: FaceVariable, *args)[source]
pyfvtool.advection.convectionUpwindTerm3D(u: FaceVariable, *args)[source]
pyfvtool.advection.convectionUpwindTermCylindrical1D(u: FaceVariable, *args)[source]
pyfvtool.advection.convectionUpwindTermCylindrical2D(u: FaceVariable, *args)[source]
pyfvtool.advection.convectionUpwindTermCylindrical3D(u: FaceVariable, *args)[source]
pyfvtool.advection.convectionUpwindTermPolar2D(u: FaceVariable, *args)[source]
pyfvtool.advection.convectionUpwindTermSpherical1D(u: FaceVariable, *args)[source]
pyfvtool.advection.convectionUpwindTermSpherical3D(u: FaceVariable, *args)[source]

Source terms

pyfvtool.source.constantSourceTerm(gamma: CellVariable)[source]

Constant source term in a PDE.

The value of this source gamma may be different in each cell, but is constant during the evolution of the PDE.

Parameters:

gamma (CellVariable) – Value of the source term

Returns:

RHS – Right hand side of the source term

Return type:

np.ndarray

Examples

>>> import pyfvtool as pf
>>> m = pf.Grid1D(10, 1.0)
>>> gamma = pf.CellVariable(m, 1.0)
>>> RHS = pf.constantSourceTerm(gamma)
pyfvtool.source.linearSourceTerm(beta: CellVariable)[source]

Linear source term in a PDE.

The linear source term takes the form beta*phi where phi is the solution variable and the factor beta a coefficient.

Parameters:

beta (CellVariable) – Multiplicative coefficient for the solution variable

Returns:

RHS – Right hand side of the source term

Return type:

np.ndarray

Examples

>>> import pyfvtool as pf
>>> m = pf.Grid1D(10, 1.0)
>>> beta = pf.CellVariable(m, 1.0)
>>> RHS = pf.linearSourceTerm(beta)
pyfvtool.source.transientTerm(phi_old: CellVariable, dt, alfa)[source]

Transient term of a PDE.

Parameters:
Returns:

  • M (csr_array) – Matrix of the transient term

  • RHS (np.ndarray) – Right hand side of the transient term

Examples

>>> import pyfvtool as pf
>>> m = pf.Grid1D(10, 1.0)
>>> phi_old = pf.CellVariable(m, 0.0)
>>> M, RHS = pf.transientTerm(phi_old, 1.0, 1.0)

Solver

pyfvtool.pdesolver.solveExplicitPDE(phi_old: CellVariable, dt: float, RHS: ndarray) CellVariable[source]

Solve the PDE using an explicit finite volume method.

Parameters:
  • phi_old (CellVariable) – Solution of the previous time step

  • dt (float) – Time step

  • RHS (np.ndarray) – Right hand side of the linear system

Returns:

phi – Solution of the PDE

Return type:

CellVariable

pyfvtool.pdesolver.solveMatrixPDE(m: MeshStructure, M: csr_array, RHS: ndarray, externalsolver=None) CellVariable[source]

Solve the PDE discretized according to the finite volume method

This ‘expert-level’ solver routine uses the M matrix and RHS right-hand side vector directly. These matrices should be constructed beforehand by combining the different M matrices and RHS vectors generated using the xxxTerm routines, including those related to the boundary conditions.

Returns a new CellVariable without changing the input (‘old’) CellVariable.

Application of boundary condition terms and updating the related ghost cells should be handled entirely by the (expert-level) user.

Parameters:
  • m (MeshStructure) – Mesh structure

  • M (csr_array) – Matrix of the linear system

  • RHS (np.ndarray) – Right hand side of the linear system

  • externalsolver (function (optional)) – If provided, use an external sparse solver via a function call having the same interface as the default built-in SuperLU solver scipy.sparse.linalg.spsolve.

Returns:

phi – Solution of the PDE (newly created CellVariable instance)

Return type:

CellVariable

pyfvtool.pdesolver.solvePDE(phi: CellVariable, eqnterms: list, externalsolver=None) CellVariable[source]

Solve a PDE using the finite volume method

This solver routine constructs the matrix equation from the FVM-discretized terms. These terms are provided by the user as a list each term being the output of a prior call to the appropriate xxxTerm() function.

The terms provided by the user should not include any terms related to the boundary conditions, as these are handled automatically by solvePDE.

The CellVariable is updated with the solution values, and a reference to this one and the same variable is returned. If the ‘old’ input CellVariable is to be conserved, it should be copied beforehand to a separate instance.

Parameters:
  • phi (CellVariable) – Solution variable subject to the boundary conditions represented by bcterm.

  • eqnterms (list) – List of matrix equation terms generated by the xxxTerm functions. The elements of the list can either be a tuple (M, RHS), a simple csr_array matrix (M) or a 1D np.ndarray (RHS). These will be added up to create both the M and RHS of the matrix equation to be solved.

  • externalsolver (function, optional) – If specified, use an external sparse solver via a function call having the same interface as the default solver scipy.sparse.linalg.spsolve. If not specified, this default built-in SuperLU solver will be used.

Raises:

TypeError – If any of the provided terms is not conform to expectations.

Returns:

The updated CellVariable.

Return type:

CellVariable

Averaging utilities

This module contains functions for averaging node values to face values.

Functions:

cell_size_array

linearMean - linear interpolation of adjacent node values to face positions

arithmeticMean - interpolate node values to face values via arithmetic averaging

geometricMean - interpolate node values to face values via geometric averaging

harmonicMean - interpolate node values to face values via harmonic averaging

upwindMean - interpolate node values to face value based on local background velocity at face position

tvdMean - NotImplemented

pyfvtool.averaging.arithmeticMean(phi: CellVariable)[source]

Interpolate a mesh-variable defined on mesh-nodes to mesh-faces by arithmetic averaging adjacent node values.

This function gets the value of the field variable phi defined over the MeshStructure and calculates the arithmetic average on the cell faces, for a uniform Mesh:

phi[face] = 1/#_adjacent_nodes * \sum(phi[adjacent nodes])
    -->
phi[face] = 0.5 * (phi[face-0.5*dx] + phi[face+0.5*dx])
    -->
phi[face] =  (cell_width[-]*phi[-] + cell_width[+]*phi[+])/(cell_width[-]+cell_width[+])
Parameters:

phi ({CellVariable object}) – Variable defined at mesh-nodes

Returns:

out – Arithmetic average of node values

Return type:

{FaceVariable object}

See also

linearMean

estimate node values to face values via linear-interpolation of node values

geometricMean

interpolate node values to face values via geometric averaging of node values

harmonicMean

interpolate node values to face values via harmonic averaging of node values

Notes

The arithmetic mean emphasizes (is sensitive to) outliers with large values.

Weighting the average by the cell length in the appropriate dimension is necessary for multi-dimensional grids, and is a necessary step in the extension to non-uniform grid spacing in 1-dimension.

Examples

>>> import pyfvtool as pf
>>> mesh = pf.Grid1D(10, 1.0)
>>> phi = pf.CellVariable(mesh, 1.0)
>>> phi_face = pf.arithmeticMean(phi)
pyfvtool.averaging.cell_size_array(m: MeshStructure)[source]
pyfvtool.averaging.geometricMean(phi: CellVariable)[source]

Interpolate a mesh-variable defined on mesh-nodes to mesh-faces by geometric averaging adjacent node values.

This function gets the value of the field variable phi defined over the MeshStructure and calculates the geometric average on the cell faces, for a uniform Mesh:

phi[face] = ( \prod(phi_nodes) )^(1/#_adjacent_nodes)
    -->
phi[face] = (phi[face-0.5*dx]*phi[face+0.5*dx])^(0.5)
    -->
phi[face] =  exp(
        ( cell_width[-] * ln(phi[-]) + cell_width[+] * ln(phi[+]) ) / ( cell_width[-] + cell_width[+] )
    )
Parameters:

phi ({CellVariable object}) – Variable defined at mesh-nodes

Returns:

out – Geometric mean of node values at face positions

Return type:

{FaceVariable object}

See also

linearMean

estimate node values to face values via linear-interpolation of node values

arithmeticMean

interpolate node values to face values via arithmetic averaging of node values

harmonicMean

interpolate node values to face values via harmonic averaging of node values

Notes

The geometric mean is equivalent to the arithmetic mean on a log-scale.
phi_face = ( prod(phi_nodes) )^(1/N_nodes)

= exp( ln( arithmeticMean(…) ) )

It is always intermediate to the harmonic or arithmetic means, and tends to the central value of th eset by ‘weight’.

Weighting the average by the cell length in the appropriate dimension is necessary for multi-dimensional grids, and is a necessary step in the extension to non-uniform grid spacing in 1-dimension.

Examples

>>>
pyfvtool.averaging.harmonicMean(phi: CellVariable)[source]

Interpolate a mesh-variable defined on mesh-nodes to mesh-faces by harmonic averaging adjacent node values.

This function gets the value of the field variable phi defined over the MeshStructure and calculates the harmonic average on the cell faces, for a uniform Mesh:

phi[face] = \sum(1.0/phi_nodes)^(-1)
    -->
phi[face] = (1/phi[face-0.5*dx] + 1/phi[face+0.5*dx])^(-1)
    -->
phi[face] = ( cell_width[-] + cell_width[+] ) / (  cell_width[-]/phi[-] + cell_width[+]/phi[+]  )
Parameters:

phi ({CellVariable object}) – Variable defined at mesh-nodes

Returns:

out – Harmonic average of node values

Return type:

{FaceVariable object}

See also

linearMean

estimate node values to face values via linear-interpolation of node values

geometricMean

interpolate node values to face values via geometric averaging of node values

arithmeticMean

interpolate node values to face values via arithmetic averaging of node values

Notes

The harmonic mean emphasizes (is sensitive to) outliers with small values.

Examples

>>>
pyfvtool.averaging.linearMean(phi: CellVariable) FaceVariable[source]

Linearly interpolate a mesh-variable defined on mesh-nodes to mesh-faces.

This function gets the value of the field variable phi defined over the MeshStructure and calculates the value on the cell faces via linear interpolation, for a uniform Mesh.

Parameters:

phi ({CellVariable object}) – Variable defined at mesh-nodes

Returns:

out – Interpolated node values at face positions

Return type:

{FaceVariable object}

See also

arithmeticMean

interpolate node values to face values via arithmetic averaging of node values

geometricMean

interpolate node values to face values via geometric averaging of node values

harmonicMean

interpolate node values to face values via harmonic averaging of node values

Notes

Examples

>>>
pyfvtool.averaging.tvdMean(phi: CellVariable, u: FaceVariable, FL)[source]

Interpolate a mesh-variable defined on mesh-nodes to mesh-faces by matching the node value carried by the flow on a grid with a Total Variation Diminishing (TVD) discretization.

This function gets the value of the field variable phi defined over the MeshStructure and determines the value on the cell faces by applying a flow-weighted average / including a linearization of the local flow (first order hold in the direction of the flow).

Parameters:
  • phi ({CellVariable object}) – Variable defined at mesh-nodes

  • u ({FaceVariable object}) – Flow-velocity at mesh-faces

Returns:

out – Interpolated value of phi at the mesh-faces

Return type:

{FaceVariable object}

See also

upwindMean

determines the direction of information flow based on the current local value of the background flow.

Notes

The TVD mean determines the direction of information flow based on a local linearization of the current local value of the background velocity. s. https://en.wikipedia.org/wiki/Total_variation_diminishing

Examples

>>>
pyfvtool.averaging.upwindMean(phi: CellVariable, u: FaceVariable)[source]

Interpolate a mesh-variable defined on mesh-nodes to mesh-faces by matching the node value carried by the flow.

This function gets the value of the field variable phi defined over the MeshStructure and determines the value on the cell faces by applying a zeroeth order hold in the direction of the flow.

phi[face] = phi[adjacent node from flow]

Parameters:
  • phi ({CellVariable object}) – Variable defined at mesh-nodes

  • u ({FaceVariable object}) – Flow-velocity at mesh-faces

Returns:

out – Zeroeth order hold in flow direction

Return type:

{FaceVariable object}

See also

linearMean

estimate node values to face values via linear-interpolation of node values

arithmeticMean

interpolate node values to face values via arithmetic averaging of node values

geometricMean

interpolate node values to face values via geometric averaging of node values

harmonicMean

interpolate node values to face values via harmonic averaging of node values

Notes

The upwind mean determines the direction of information flow based on the current local value of the background flow.

Examples

>>>

Visualization

pyfvtool.visualization.visualizeCells(phi: CellVariable, vmin=None, vmax=None, cmap='viridis', shading='gouraud')[source]

Visualize the cell variable.

Parameters:
  • phi (CellVariable) – Cell variable to be visualized

  • vmin (float) – Minimum value of the colormap

  • vmax (float) – Maximum value of the colormap

  • cmap (str) – Colormap

  • shading (str) – Shading method

Examples

>>> import pyfvtool as pf
>>> m = pf.Grid1D(10, 1.0)
>>> phi = pf.CellVariable(m, 1.0)
>>> pf.visualizeCells(phi)