discretize.TreeMesh

class discretize.TreeMesh(*args, **kwargs)[source]

Bases: discretize.tree_ext._TreeMesh, discretize.base.base_tensor_mesh.BaseTensorMesh, discretize.InnerProducts.InnerProducts, discretize.MeshIO.TreeMeshIO

TreeMesh is a class for adaptive QuadTree (2D) and OcTree (3D) meshes.

Required Properties:

  • axis_u (Vector3): Vector orientation of u-direction. For more details see the docs for the rotation_matrix property., a 3D Vector of <class ‘float’> with shape (3), Default: X

  • axis_v (Vector3): Vector orientation of v-direction. For more details see the docs for the rotation_matrix property., a 3D Vector of <class ‘float’> with shape (3), Default: Y

  • axis_w (Vector3): Vector orientation of w-direction. For more details see the docs for the rotation_matrix property., a 3D Vector of <class ‘float’> with shape (3), Default: Z

  • h (a list of Array): h is a list containing the cell widths of the tensor mesh in each dimension., a list (each item is a list or numpy array of <class ‘float’> with shape (*)) with length between 0 and 3

  • reference_system (String): The type of coordinate reference frame. Can take on the values cartesian, cylindrical, or spherical. Abbreviations of these are allowed., a unicode string, Default: cartesian

  • x0 (Array): origin of the mesh (dim, ), a list or numpy array of <class ‘float’>, <class ‘int’> with shape (*)

Attributes
area

Returns a numpy array of length nF with the area (length in 2D) of all faces ordered by x, then y, then z.

aveCC2F

Construct the averaging operator on cell centers to cell faces.

aveCC2Fx

Construct the averaging operator on cell centers to cell x-faces.

aveCC2Fy

Construct the averaging operator on cell centers to cell y-faces.

aveCC2Fz

Construct the averaging operator on cell centers to cell z-faces.

aveCCV2F

Construct the averaging operator on cell centers to cell faces.

aveE2CC

Construct the averaging operator on cell edges to cell centers.

aveE2CCV

Construct the averaging operator on cell edges to cell centers.

aveEx2CC

Construct the averaging operator on cell edges in the x direction to cell centers.

aveEy2CC

Construct the averaging operator on cell edges in the y direction to cell centers.

aveEz2CC

Construct the averaging operator on cell edges in the z direction to cell centers.

aveF2CC

Construct the averaging operator on cell faces to cell centers.

aveF2CCV

Construct the averaging operator on cell faces to cell centers.

aveFx2CC

Construct the averaging operator on cell faces in the x direction to cell centers.

aveFy2CC

Construct the averaging operator on cell faces in the y direction to cell centers.

aveFz2CC

Construct the averaging operator on cell faces in the z direction to cell centers.

aveN2CC

Construct the averaging operator on cell nodes to cell centers.

aveN2E

Construct the averaging operator on cell nodes to cell edges, keeping each dimension separate.

aveN2Ex

Averaging operator on cell nodes to x-edges

aveN2Ey

Averaging operator on cell nodes to y-edges

aveN2Ez

Averaging operator on cell nodes to z-edges

aveN2F

Construct the averaging operator on cell nodes to cell edges, keeping each dimension separate.

aveN2Fx

Averaging operator on cell nodes to x-faces

aveN2Fy

Averaging operator on cell nodes to y-faces

aveN2Fz

Averaging operator on cell nodes to z-faces

axis_u

axis_u (Vector3): Vector orientation of u-direction. For more details see the docs for the rotation_matrix property., a 3D Vector of <class ‘float’> with shape (3), Default: X

axis_v

axis_v (Vector3): Vector orientation of v-direction. For more details see the docs for the rotation_matrix property., a 3D Vector of <class ‘float’> with shape (3), Default: Y

axis_w

axis_w (Vector3): Vector orientation of w-direction. For more details see the docs for the rotation_matrix property., a 3D Vector of <class ‘float’> with shape (3), Default: Z

cellBoundaryInd

Returns a tuple of arrays of indexes for boundary cells in each direction

cellGrad

Cell centered Gradient operator built off of the faceDiv operator.

cellGradStencil
cellGradx

Cell centered Gradient operator in x-direction (Gradx)

cellGrady

Cell centered Gradient operator in y-direction (Grady)

cellGradz

Cell centered Gradient operator in z-direction (Gradz)

dim

The dimension of the mesh (1, 2, or 3).

edge

Returns a numpy array of length nE with the length of all edges ordered by x, then y, then z.

edgeCurl

Construct the 3D curl operator.

faceBoundaryInd

Returns a tuple of arrays of indexes for boundary faces in each direction

faceDiv

Construct divergence operator (face-stg to cell-centres).

faceDivx
faceDivy
faceDivz
fill

How filled is the mesh compared to a TensorMesh? As a fraction: [0, 1].

gridCC

Returns a numpy arrayof shape (nC, dim) with the center locations of all cells in order.

gridEx

Returns a numpy array of shape (nEx, dim) with the centers of all non-hanging edges along the first dimension in order.

gridEy

Returns a numpy array of shape (nEy, dim) with the centers of all non-hanging edges along the second dimension in order.

gridEz

Returns a numpy array of shape (nEz, dim) with the centers of all non-hanging edges along the third dimension in order.

gridFx

Returns a numpy array of shape (nFx, dim) with the centers of all non-hanging faces along the first dimension in order.

gridFy

Returns a numpy array of shape (nFy, dim) with the centers of all non-hanging faces along the second dimension in order.

gridFz

Returns a numpy array of shape (nFz, dim) with the centers of all non-hanging faces along the third dimension in order.

gridN

Returns a numpy array of shape (nN, dim) with the locations of all non-hanging nodes in order.

gridhEx

Returns a numpy array of shape (nhEx, dim) with the centers of all hanging edges along the first dimension in order.

gridhEy

Returns a numpy array of shape (nhEy, dim) with the centers of all hanging edges along the second dimension in order.

gridhEz

Returns a numpy array of shape (nhEz, dim) with the centers of all hanging edges along the third dimension in order.

gridhFx

Returns a numpy array of shape (nhFx, dim) with the centers of all hanging faces along the first dimension in order.

gridhFy

Returns a numpy array of shape (nhFy, dim) with the centers of all hanging faces along the second dimension in order.

gridhFz

Returns a numpy array of shape (nhFz, dim) with the centers of all hanging faces along the third dimension in order.

gridhN

Returns a numpy array of shape (nN, dim) with the locations of all hanging nodes in order.

h

h (a list of Array): h is a list containing the cell widths of the tensor mesh in each dimension., a list (each item is a list or numpy array of <class ‘float’> with shape (*)) with length between 0 and 3

h_gridded

Returns an (nC, dim) numpy array with the widths of all cells in order

hx

Width of cells in the x direction

hy

Width of cells in the y direction

hz

Width of cells in the z direction

maxLevel

The maximum level used, which may be less than levels or max_level.

max_level

The maximum possible level for a cell on this mesh

nC

Number of cells

nE

Total number of non-hanging edges amongst all dimensions

nEx

Number of non-hanging edges oriented along the first dimension

nEy

Number of non-hanging edges oriented along the second dimension

nEz

Number of non-hanging edges oriented along the third dimension

nF

Total number of non-hanging faces amongst all dimensions

nFx

Number of non-hanging faces oriented along the first dimension

nFy

Number of non-hanging faces oriented along the second dimension

nFz

Number of non-hanging faces oriented along the third dimension

nN

Number of non-hanging nodes

nhE

Total number of hanging edges amongst all dimensions

nhEx

Number of hanging edges oriented along the first dimension

nhEy

Number of hanging edges oriented along the second dimension

nhEz

Number of hanging edges oriented along the third dimension

nhF

Total number of hanging faces amongst all dimensions

nhFx

Number of hanging faces oriented along the first dimension

nhFy

Number of hanging faces oriented along the second dimension

nhFz

Number of hanging faces oriented along the third dimension

nhN

Number of hanging nodes

nodalGrad

Construct gradient operator (nodes to edges).

nodalLaplacian
normals

Face Normals

ntE

Total number of non-hanging and hanging edges amongst all dimensions

ntEx

Number of non-hanging and hanging edges oriented along the first dimension

ntEy

Number of non-hanging and hanging edges oriented along the second dimension

ntEz

Number of non-hanging and hanging edges oriented along the third dimension

ntF

Total number of hanging and non-hanging faces amongst all dimensions

ntFx

Number of non-hanging and hanging faces oriented along the first dimension

ntFy

Number of non-hanging and hanging faces oriented along the second dimension

ntFz

Number of non-hanging and hanging faces oriented along the third dimension

ntN

Number of non-hanging and hanging nodes

permuteCC

Permutation matrix re-ordering of cells sorted by x, then y, then z

permuteE

Permutation matrix re-ordering of edges sorted by x, then y, then z

permuteF

Permutation matrix re-ordering of faces sorted by x, then y, then z

reference_is_rotated

True if the axes are rotated from the traditional <X,Y,Z> system

reference_system

reference_system (String): The type of coordinate reference frame. Can take on the values cartesian, cylindrical, or spherical. Abbreviations of these are allowed., a unicode string, Default: cartesian

rotation_matrix

Builds a rotation matrix to transform coordinates from their coordinate system into a conventional cartesian system.

tangents

Edge Tangents

vectorCCx

Cell-centered grid vector (1D) in the x direction.

vectorCCy

Cell-centered grid vector (1D) in the y direction.

vectorCCz

Cell-centered grid vector (1D) in the z direction.

vectorNx

Nodal grid vector (1D) in the x direction.

vectorNy

Nodal grid vector (1D) in the y direction.

vectorNz

Nodal grid vector (1D) in the z direction.

vnE

Total number of edges in each direction

vnF

Total number of faces in each direction

vntE

Total number of hanging and non-hanging edges in a [nx,ny,nz] form

vntF

Total number of hanging and non-hanging faces in a [nx,ny,nz] form

vol

Returns a numpy array of length nC with the volumes (areas in 2D) of all cells in order.

x0

x0 (Array): origin of the mesh (dim, ), a list or numpy array of <class ‘float’>, <class ‘int’> with shape (*)

Methods

copy

Generic (shallow and deep) copying operations.

plotSlice

Attributes

TreeMesh.area

Returns a numpy array of length nF with the area (length in 2D) of all faces ordered by x, then y, then z.

TreeMesh.aveCC2F

Construct the averaging operator on cell centers to cell faces.

TreeMesh.aveCC2Fx

Construct the averaging operator on cell centers to cell x-faces.

TreeMesh.aveCC2Fy

Construct the averaging operator on cell centers to cell y-faces.

TreeMesh.aveCC2Fz

Construct the averaging operator on cell centers to cell z-faces.

TreeMesh.aveCCV2F

Construct the averaging operator on cell centers to cell faces.

TreeMesh.aveE2CC

Construct the averaging operator on cell edges to cell centers.

TreeMesh.aveE2CCV

Construct the averaging operator on cell edges to cell centers.

TreeMesh.aveEx2CC

Construct the averaging operator on cell edges in the x direction to cell centers.

TreeMesh.aveEy2CC

Construct the averaging operator on cell edges in the y direction to cell centers.

TreeMesh.aveEz2CC

Construct the averaging operator on cell edges in the z direction to cell centers.

TreeMesh.aveF2CC

Construct the averaging operator on cell faces to cell centers.

TreeMesh.aveF2CCV

Construct the averaging operator on cell faces to cell centers.

TreeMesh.aveFx2CC

Construct the averaging operator on cell faces in the x direction to cell centers.

TreeMesh.aveFy2CC

Construct the averaging operator on cell faces in the y direction to cell centers.

TreeMesh.aveFz2CC

Construct the averaging operator on cell faces in the z direction to cell centers.

TreeMesh.aveN2CC

Construct the averaging operator on cell nodes to cell centers.

TreeMesh.aveN2E

Construct the averaging operator on cell nodes to cell edges, keeping each dimension separate.

TreeMesh.aveN2Ex

Averaging operator on cell nodes to x-edges

TreeMesh.aveN2Ey

Averaging operator on cell nodes to y-edges

TreeMesh.aveN2Ez

Averaging operator on cell nodes to z-edges

TreeMesh.aveN2F

Construct the averaging operator on cell nodes to cell edges, keeping each dimension separate.

TreeMesh.aveN2Fx

Averaging operator on cell nodes to x-faces

TreeMesh.aveN2Fy

Averaging operator on cell nodes to y-faces

TreeMesh.aveN2Fz

Averaging operator on cell nodes to z-faces

TreeMesh.axis_u

axis_u (Vector3): Vector orientation of u-direction. For more details see the docs for the rotation_matrix property., a 3D Vector of <class ‘float’> with shape (3), Default: X

TreeMesh.axis_v

axis_v (Vector3): Vector orientation of v-direction. For more details see the docs for the rotation_matrix property., a 3D Vector of <class ‘float’> with shape (3), Default: Y

TreeMesh.axis_w

axis_w (Vector3): Vector orientation of w-direction. For more details see the docs for the rotation_matrix property., a 3D Vector of <class ‘float’> with shape (3), Default: Z

TreeMesh.cellBoundaryInd

Returns a tuple of arrays of indexes for boundary cells in each direction xdown, xup, ydown, yup, zdown, zup

TreeMesh.cellGrad

Cell centered Gradient operator built off of the faceDiv operator. Grad = - (Mf)^{-1} * Div * diag (volume)

TreeMesh.cellGradStencil
TreeMesh.cellGradx

Cell centered Gradient operator in x-direction (Gradx) Grad = sp.vstack((Gradx, Grady, Gradz))

TreeMesh.cellGrady

Cell centered Gradient operator in y-direction (Grady) Grad = sp.vstack((Gradx, Grady, Gradz))

TreeMesh.cellGradz

Cell centered Gradient operator in z-direction (Gradz) Grad = sp.vstack((Gradx, Grady, Gradz))

TreeMesh.dim

The dimension of the mesh (1, 2, or 3).

Returns
int

dimension of the mesh

TreeMesh.edge

Returns a numpy array of length nE with the length of all edges ordered by x, then y, then z.

TreeMesh.edgeCurl

Construct the 3D curl operator.

TreeMesh.faceBoundaryInd

Returns a tuple of arrays of indexes for boundary faces in each direction xdown, xup, ydown, yup, zdown, zup

TreeMesh.faceDiv

Construct divergence operator (face-stg to cell-centres).

TreeMesh.faceDivx
TreeMesh.faceDivy
TreeMesh.faceDivz
TreeMesh.fill

How filled is the mesh compared to a TensorMesh? As a fraction: [0, 1].

TreeMesh.gridCC

Returns a numpy arrayof shape (nC, dim) with the center locations of all cells in order.

TreeMesh.gridEx

Returns a numpy array of shape (nEx, dim) with the centers of all non-hanging edges along the first dimension in order.

TreeMesh.gridEy

Returns a numpy array of shape (nEy, dim) with the centers of all non-hanging edges along the second dimension in order.

TreeMesh.gridEz

Returns a numpy array of shape (nEz, dim) with the centers of all non-hanging edges along the third dimension in order.

TreeMesh.gridFx

Returns a numpy array of shape (nFx, dim) with the centers of all non-hanging faces along the first dimension in order.

TreeMesh.gridFy

Returns a numpy array of shape (nFy, dim) with the centers of all non-hanging faces along the second dimension in order.

TreeMesh.gridFz

Returns a numpy array of shape (nFz, dim) with the centers of all non-hanging faces along the third dimension in order.

TreeMesh.gridN

Returns a numpy array of shape (nN, dim) with the locations of all non-hanging nodes in order.

TreeMesh.gridhEx

Returns a numpy array of shape (nhEx, dim) with the centers of all hanging edges along the first dimension in order.

TreeMesh.gridhEy

Returns a numpy array of shape (nhEy, dim) with the centers of all hanging edges along the second dimension in order.

TreeMesh.gridhEz

Returns a numpy array of shape (nhEz, dim) with the centers of all hanging edges along the third dimension in order.

TreeMesh.gridhFx

Returns a numpy array of shape (nhFx, dim) with the centers of all hanging faces along the first dimension in order.

TreeMesh.gridhFy

Returns a numpy array of shape (nhFy, dim) with the centers of all hanging faces along the second dimension in order.

TreeMesh.gridhFz

Returns a numpy array of shape (nhFz, dim) with the centers of all hanging faces along the third dimension in order.

TreeMesh.gridhN

Returns a numpy array of shape (nN, dim) with the locations of all hanging nodes in order.

TreeMesh.h

h (a list of Array): h is a list containing the cell widths of the tensor mesh in each dimension., a list (each item is a list or numpy array of <class ‘float’> with shape (*)) with length between 0 and 3

TreeMesh.h_gridded

Returns an (nC, dim) numpy array with the widths of all cells in order

TreeMesh.hx

Width of cells in the x direction

TreeMesh.hy

Width of cells in the y direction

TreeMesh.hz

Width of cells in the z direction

TreeMesh.maxLevel

The maximum level used, which may be less than levels or max_level.

TreeMesh.max_level

The maximum possible level for a cell on this mesh

TreeMesh.nC

Number of cells

TreeMesh.nE

Total number of non-hanging edges amongst all dimensions

TreeMesh.nEx

Number of non-hanging edges oriented along the first dimension

TreeMesh.nEy

Number of non-hanging edges oriented along the second dimension

TreeMesh.nEz

Number of non-hanging edges oriented along the third dimension

TreeMesh.nF

Total number of non-hanging faces amongst all dimensions

TreeMesh.nFx

Number of non-hanging faces oriented along the first dimension

TreeMesh.nFy

Number of non-hanging faces oriented along the second dimension

TreeMesh.nFz

Number of non-hanging faces oriented along the third dimension

TreeMesh.nN

Number of non-hanging nodes

TreeMesh.nhE

Total number of hanging edges amongst all dimensions

TreeMesh.nhEx

Number of hanging edges oriented along the first dimension

TreeMesh.nhEy

Number of hanging edges oriented along the second dimension

TreeMesh.nhEz

Number of hanging edges oriented along the third dimension

TreeMesh.nhF

Total number of hanging faces amongst all dimensions

TreeMesh.nhFx

Number of hanging faces oriented along the first dimension

TreeMesh.nhFy

Number of hanging faces oriented along the second dimension

TreeMesh.nhFz

Number of hanging faces oriented along the third dimension

TreeMesh.nhN

Number of hanging nodes

TreeMesh.nodalGrad

Construct gradient operator (nodes to edges).

TreeMesh.nodalLaplacian
TreeMesh.normals

Face Normals

Return type

numpy.ndarray

Returns

normals, (sum(nF), dim)

TreeMesh.ntE

Total number of non-hanging and hanging edges amongst all dimensions

TreeMesh.ntEx

Number of non-hanging and hanging edges oriented along the first dimension

TreeMesh.ntEy

Number of non-hanging and hanging edges oriented along the second dimension

TreeMesh.ntEz

Number of non-hanging and hanging edges oriented along the third dimension

TreeMesh.ntF

Total number of hanging and non-hanging faces amongst all dimensions

TreeMesh.ntFx

Number of non-hanging and hanging faces oriented along the first dimension

TreeMesh.ntFy

Number of non-hanging and hanging faces oriented along the second dimension

TreeMesh.ntFz

Number of non-hanging and hanging faces oriented along the third dimension

TreeMesh.ntN

Number of non-hanging and hanging nodes

TreeMesh.permuteCC

Permutation matrix re-ordering of cells sorted by x, then y, then z

TreeMesh.permuteE

Permutation matrix re-ordering of edges sorted by x, then y, then z

TreeMesh.permuteF

Permutation matrix re-ordering of faces sorted by x, then y, then z

TreeMesh.reference_is_rotated

True if the axes are rotated from the traditional <X,Y,Z> system with vectors of \((1,0,0)\), \((0,1,0)\), and \((0,0,1)\)

TreeMesh.reference_system

reference_system (String): The type of coordinate reference frame. Can take on the values cartesian, cylindrical, or spherical. Abbreviations of these are allowed., a unicode string, Default: cartesian

TreeMesh.rotation_matrix

Builds a rotation matrix to transform coordinates from their coordinate system into a conventional cartesian system. This is built off of the three axis_u, axis_v, and axis_w properties; these mapping coordinates use the letters U, V, and W (the three letters preceding X, Y, and Z in the alphabet) to define the projection of the X, Y, and Z durections. These UVW vectors describe the placement and transformation of the mesh’s coordinate sytem assuming at most 3 directions.

Why would you want to use these UVW mapping vectors the this rotation_matrix property? They allow us to define the realationship between local and global coordinate systems and provide a tool for switching between the two while still maintaing the connectivity of the mesh’s cells. For a visual example of this, please see the figure in the docs for the InterfaceVTK.

TreeMesh.tangents

Edge Tangents

Return type

numpy.ndarray

Returns

normals, (sum(nE), dim)

TreeMesh.vectorCCx

Cell-centered grid vector (1D) in the x direction.

TreeMesh.vectorCCy

Cell-centered grid vector (1D) in the y direction.

TreeMesh.vectorCCz

Cell-centered grid vector (1D) in the z direction.

TreeMesh.vectorNx

Nodal grid vector (1D) in the x direction.

TreeMesh.vectorNy

Nodal grid vector (1D) in the y direction.

TreeMesh.vectorNz

Nodal grid vector (1D) in the z direction.

TreeMesh.vnE

Total number of edges in each direction

Returns
vnEnumpy.ndarray = [nEx, nEy, nEz], (dim, )
import discretize
import numpy as np
M = discretize.TensorMesh([np.ones(n) for n in [2,3]])
M.plotGrid(edges=True, show_it=True)

(Source code, png, pdf)

../../_images/discretize-TreeMesh-1.png
TreeMesh.vnF

Total number of faces in each direction

Return type

numpy.ndarray

Returns

[nFx, nFy, nFz], (dim, )

import discretize
import numpy as np
M = discretize.TensorMesh([np.ones(n) for n in [2,3]])
M.plotGrid(faces=True, show_it=True)

(Source code, png, pdf)

../../_images/discretize-TreeMesh-2.png
TreeMesh.vntE

Total number of hanging and non-hanging edges in a [nx,ny,nz] form

TreeMesh.vntF

Total number of hanging and non-hanging faces in a [nx,ny,nz] form

TreeMesh.vol

Returns a numpy array of length nC with the volumes (areas in 2D) of all cells in order.

TreeMesh.x0

x0 (Array): origin of the mesh (dim, ), a list or numpy array of <class ‘float’>, <class ‘int’> with shape (*)

Methods

TreeMesh.cell_levels_by_index(indices)[source]

Fast function to return a list of levels for the given cell indices

Parameters
index: array_like of length (N)

Cell indexes to query

Returns
numpy.array of length (N)

Levels for the cells.

TreeMesh.copy()

Make a copy of the current mesh

classmethod TreeMesh.deserialize(serial)[source]

Creates HasProperties instance from serialized dictionary

This uses the Property deserializers to deserialize all JSON-compatible dictionary values into their corresponding Property values on a new instance of a HasProperties class. Extra keys in the dictionary that do not correspond to Properties will be ignored.

Parameters:

  • value - Dictionary to deserialize new instance from.

  • trusted - If True (and if the input dictionary has '__class__' keyword and this class is in the registry), the new HasProperties class will come from the dictionary. If False (the default), only the HasProperties class this method is called on will be constructed.

  • strict - Requires '__class__', if present on the input dictionary, to match the deserialized instance’s class. Also disallows unused properties in the input dictionary. Default is False.

  • assert_valid - Require deserialized instance to be valid. Default is False.

  • Any other keyword arguments will be passed through to the Property deserializers.

TreeMesh.equal(other)

Determine if two HasProperties instances are equivalent

Equivalence is determined by checking if all Property values on two instances are equal, using Property.equal.

TreeMesh.finalize(self)

Finalize the TreeMesh Called after finished cronstruction of the mesh. Can only be called once. After finalize is called, all other attributes and functions are valid.

static TreeMesh.from_omf(element)

Convert an OMF element to it’s proper discretize type. Automatically determines the output type. Returns both the mesh and a dictionary of model arrays.

TreeMesh.getEdgeInnerProduct(prop=None, invProp=False, invMat=False, doFast=True)

Generate the edge inner product matrix

Parameters
propnumpy.ndarray

material property (tensor properties are possible) at each cell center (nC, (1, 3, or 6))

invPropbool

inverts the material property

invMatbool

inverts the matrix

doFastbool

do a faster implementation if available.

Returns
scipy.sparse.csr_matrix

M, the inner product matrix (nE, nE)

TreeMesh.getEdgeInnerProductDeriv(prop, doFast=True, invProp=False, invMat=False)
Parameters
propnumpy.ndarray

material property (tensor properties are possible) at each cell center (nC, (1, 3, or 6))

doFastbool

do a faster implementation if available.

invPropbool

inverts the material property

invMatbool

inverts the matrix

Returns
scipy.sparse.csr_matrix

dMdm, the derivative of the inner product matrix (nE, nC*nA)

TreeMesh.getFaceInnerProduct(prop=None, invProp=False, invMat=False, doFast=True)

Generate the face inner product matrix

Parameters
propnumpy.ndarray

material property (tensor properties are possible) at each cell center (nC, (1, 3, or 6))

invPropbool

inverts the material property

invMatbool

inverts the matrix

doFastbool

do a faster implementation if available.

Returns
scipy.sparse.csr_matrix

M, the inner product matrix (nF, nF)

TreeMesh.getFaceInnerProductDeriv(prop, doFast=True, invProp=False, invMat=False)
Parameters
propnumpy.ndarray

material property (tensor properties are possible) at each cell center (nC, (1, 3, or 6))

doFast :

bool do a faster implementation if available.

invPropbool

inverts the material property

invMatbool

inverts the matrix

Returns
scipy.sparse.csr_matrix

dMdmu(u), the derivative of the inner product matrix for a certain u

TreeMesh.getInterpolationMat(locs, locType, zerosOutside=False)[source]

Produces interpolation matrix

Parameters
locnumpy.ndarray

Location of points to interpolate to

locType: str

What to interpolate

locType can be:

'Ex'    -> x-component of field defined on edges
'Ey'    -> y-component of field defined on edges
'Ez'    -> z-component of field defined on edges
'Fx'    -> x-component of field defined on faces
'Fy'    -> y-component of field defined on faces
'Fz'    -> z-component of field defined on faces
'N'     -> scalar field defined on nodes
'CC'    -> scalar field defined on cell centers
Returns
scipy.sparse.csr_matrix

M, the interpolation matrix

TreeMesh.getTensor(key)

Returns a tensor list.

Parameters
keystr

Which tensor (see below)

key can be:

'CC'    -> scalar field defined on cell centers
'N'     -> scalar field defined on nodes
'Fx'    -> x-component of field defined on faces
'Fy'    -> y-component of field defined on faces
'Fz'    -> z-component of field defined on faces
'Ex'    -> x-component of field defined on edges
'Ey'    -> y-component of field defined on edges
'Ez'    -> z-component of field defined on edges
Returns
list

list of the tensors that make up the mesh.

TreeMesh.get_boundary_cells(self, active_ind=None, direction='zu')

Returns the indices of boundary cells in a given direction given an active index array.

Parameters
active_indarray_like of bool, optional

If not None, then this must show which cells are active

direction: str, optional

must be one of (‘zu’, ‘zd’, ‘xu’, ‘xd’, ‘yu’, ‘yd’)

Returns
numpy.array

Array of indices for the boundary cells in the requested direction

TreeMesh.get_cells_along_line(self, x0, x1)

Finds the cells along a line segment defined by two points

Parameters
x0,x1array_like of length (dim)

Begining and ending point of the line segment.

Returns
list of ints

Indexes for cells that contain the a line defined by the two input points, ordered in the direction of the line.

TreeMesh.get_overlapping_cells(self, rectangle)
TreeMesh.insert_cells(self, points, levels, finalize=True)

Insert cells into the TreeMesh that contain given points

Insert cell(s) into the TreeMesh that contain the given point(s) at the assigned level(s).

Parameters
pointsarray_like with shape (N, dim)
levelsarray_like of integers with shape (N)
finalizebool, optional

Whether to finalize after inserting point(s)

Examples

>>> from discretize import TreeMesh
>>> mesh = TreeMesh([32,32])
>>> mesh.insert_cells([0.5, 0.5], mesh.max_level)
>>> print(mesh)
---- QuadTreeMesh ----
 x0: 0.00
 y0: 0.00
 hx: 32*0.03,
 hy: 32*0.03,
nC: 40
Fill: 3.91%
TreeMesh.isInside(pts, locType='N')

Determines if a set of points are inside a mesh.

Parameters

pts (numpy.ndarray) – Location of points to test

Return type

numpy.ndarray

Returns

inside, numpy array of booleans

TreeMesh.number(self)

Number the cells, nodes, faces, and edges of the TreeMesh

TreeMesh.plotGrid(self, ax=None, nodes=False, faces=False, centers=False, edges=False, lines=True, cell_line=False, faces_x=False, faces_y=False, faces_z=False, edges_x=False, edges_y=False, edges_z=False, show_it=False, **kwargs)

Plot the nodel, cell-centered, and staggered grids for 2 and 3 dimensions Plots the mesh grid in either 2D or 3D of the TreeMesh

Parameters
axmatplotlib.axes.Axes or None, optional

The axes handle to plot on

nodesbool, optional

Plot the nodal points

facesbool, optional

Plot the center points of the faces

centersbool, optional

Plot the center points of the cells

edgesbool, optional

Plot the center points of the edges

linesbool, optional

Plot the lines connecting the nodes

cell_linebool, optional

Plot the line through the cell centers in order

faces_x, faces_y, faces_zbool, optional

Plot the center points of the x, y, or z faces

edges_x, edges_y, edges_zbool, optional

Plot the center points of the x, y, or z edges

show_itbool, optional

whether to call plt.show() within the codes

colorColor or str, optional

if lines=True, the color of the lines, defaults to first color.

linewidthfloat, optional

if lines=True, the linewidth for the lines.

Returns
matplotlib.axes.Axes

Axes handle for the plot

TreeMesh.plotImage(self, v, v_type='CC', grid=False, view='real', ax=None, clim=None, show_it=False, pcolor_opts=None, grid_opts=None, range_x=None, range_y=None, **kwargs)

Plots an image of values defined on the TreeMesh If 3D, this function plots a default slice of the TreeMesh

Parameters
varray_like

Array containing the values to plot

v_typestr, optional

type of value in v one of ‘CC’, ‘N’, ‘Fx’, ‘Fy’, ‘Fz’, ‘Ex’, ‘Ey’, or ‘Ez’

gridbool, optional

plot the grid lines

view[‘real’, ‘imag’, ‘abs’], optional

The values to plot from v

axmatplotlib.axes.Axes or None, optional

The axes handle

climarray_like of length 2, or None, optional

A pair of [min, max] for the Colorbar

pcolor_optsdict, or None

options to be passed on to pcolormesh

gridOptdict, or None

options for the plotting the grid

range_x, range_y: array_like, optional

pairs of [min, max] values for the x and y ranges

TreeMesh.plotSlice(v, v_type='CC', normal='Z', ind=None, grid=False, view='real', ax=None, clim=None, show_it=False, pcolor_opts=None, stream_opts=None, grid_opts=None, range_x=None, range_y=None, **kwargs)[source]
TreeMesh.point2index(locs)[source]

Finds cells that contain the given points. Returns an array of index values of the cells that contain the given points

Parameters
locs: array_like of shape (N, dim)

points to search for the location of

Returns
numpy.array of integers of length(N)

Cell indices that contain the points

TreeMesh.projectEdgeVector(eV)

Given a vector, eV, in cartesian coordinates, this will project it onto the mesh using the tangents

Parameters

eV (numpy.ndarray) – edge vector with shape (nE, dim)

Return type

numpy.ndarray

Returns

projected edge vector, (nE, )

TreeMesh.projectFaceVector(fV)

Given a vector, fV, in cartesian coordinates, this will project it onto the mesh using the normals

Parameters

fV (numpy.ndarray) – face vector with shape (nF, dim)

Return type

numpy.ndarray

Returns

projected face vector, (nF, )

TreeMesh.readModelUBC(fileName)

Read UBC OcTree model and get vector :param string fileName: path to the UBC GIF model file to read :rtype: numpy.ndarray :return: OcTree model

classmethod TreeMesh.readUBC(meshFile)

Read UBC 3D OcTree mesh file Input: :param str meshFile: path to the UBC GIF OcTree mesh file to read :rtype: discretize.TreeMesh :return: The octree mesh

TreeMesh.refine(self, function, finalize=True)

Refine a TreeMesh using a user supplied function.

Refines the TreeMesh using a function that is recursively called on each cell of the mesh. It must accept an object of type discretize.TreeMesh.Cell and return an integer like object defining the desired level. The function can also simply be an integer, which will then cause all cells to be at least that level.

Parameters
functioncallable | int

a function describing the desired level, or an integer to refine all cells to at least that level.

finalizebool, optional

Whether to finalize the mesh

See also

discretize.TreeMesh.TreeCell

a description of the TreeCell object

Examples

>>> from discretize import TreeMesh
>>> mesh = TreeMesh([32,32])
>>> def func(cell):
>>>     r = np.linalg.norm(cell.center-0.5)
>>>     return mesh.max_level if r<0.2 else mesh.max_level-1
>>> mesh.refine(func)
>>> mesh
---- QuadTreeMesh ----
 x0: 0.00
 y0: 0.00
 hx: 32*0.03,
 hy: 32*0.03,
nC: 352
Fill: 34.38%
TreeMesh.save(filename='mesh.json', verbose=False)

Save the mesh to json :param str file: filename for saving the casing properties :param str directory: working directory for saving the file

TreeMesh.serialize()[source]

Serializes a HasProperties instance to dictionary

This uses the Property serializers to serialize all Property values to a JSON-compatible dictionary. Properties that are undefined are not included. If the HasProperties instance contains a reference to itself, a properties.SelfReferenceError will be raised.

Parameters:

  • include_class - If True (the default), the name of the class will also be saved to the serialized dictionary under key '__class__'

  • save_dynamic - If True, dynamic properties are written to the serialized dict (default: False).

  • Any other keyword arguments will be passed through to the Property serializers.

TreeMesh.toVTK(models=None)

Convert this mesh object to it’s proper VTK or pyvista data object with the given model dictionary as the cell data of that dataset.

Parameters
modelsdict(numpy.ndarray)

Name(‘s) and array(‘s). Match number of cells

TreeMesh.to_omf(models=None)

Convert this mesh object to it’s proper omf data object with the given model dictionary as the cell data of that dataset.

Parameters
modelsdict(numpy.ndarray)

Name(‘s) and array(‘s). Match number of cells

TreeMesh.to_vtk(models=None)

Convert this mesh object to it’s proper VTK or pyvista data object with the given model dictionary as the cell data of that dataset.

Parameters
modelsdict(numpy.ndarray)

Name(‘s) and array(‘s). Match number of cells

TreeMesh.validate()

Call all registered class validator methods

These are all methods decorated with @properties.validator. Validator methods are expected to raise a ValidationError if they fail.

TreeMesh.writeModelUBC(fileName, model, directory='')

Writes a model associated with a TreeMesh to a UBC-GIF format model file.

Input: :param str fileName: File to write to or just its name if directory is specified :param str directory: directory where the UBC GIF file lives :param numpy.ndarray model: The model

TreeMesh.writeUBC(fileName, models=None, directory='')

Write UBC ocTree mesh and model files from a octree mesh and model. :param string fileName: File to write to :param dict models: Models in a dict, where each key is the filename :param str directory: directory where to save model(s)

TreeMesh.writeVTK(filename, models=None, directory='')

Makes and saves a VTK object from this mesh and given models

Parameters
filenamestr

path to the output vtk file or just its name if directory is specified

modelsdict

dictionary of numpy.array - Name(‘s) and array(‘s). Match number of cells

directorystr

directory where the UBC GIF file lives

TreeMesh.write_vtk(filename, models=None, directory='')

Makes and saves a VTK object from this mesh and given models

Parameters
filenamestr

path to the output vtk file or just its name if directory is specified

modelsdict

dictionary of numpy.array - Name(‘s) and array(‘s). Match number of cells

directorystr

directory where the UBC GIF file lives